First tutorial on GW¶
The quasi-particle band structure of Silicon in the GW approximation.¶
This tutorial aims at showing how to calculate self-energy corrections to the DFT Kohn-Sham (KS) eigenvalues in the GW approximation.
A brief description of the formalism and of the equations implemented in the code can be found in the GW_notes. The different formulas of the GW formalism have been written in a pdf document by Valerio Olevano who also wrote the first version of this tutorial. For a much more consistent discussion of the theoretical aspects of the GW method we refer the reader to the review article Quasiparticle calculations in solids by W.G Aulbur et al also available here.
It is suggested to acknowledge the efforts of developers of the GW part of ABINIT, by citing the 2005 ABINIT publication.
The user should be familiarized with the four basic tutorials of ABINIT, see the tutorial home page After this first tutorial on GW, you should read the second GW tutorial.
This tutorial should take about 2 hours.
Note
Supposing you made your own install of ABINIT, the input files to run the examples are in the ~abinit/tests/ directory where ~abinit is the absolute path of the abinit top-level directory. If you have NOT made your own install, ask your system administrator where to find the package, especially the executable and test files.
To execute the tutorials, create a working directory (Work*
) and
copy there the input files and the files file of the lesson. This will be explicitly mentioned in the first lessons,
that will tell you more about the files file (see also section 1.1).
The files file ending with _x (e.g. tbase1_x.files) must be edited every time you start to use a new input file.
Most of the tutorials do not rely on parallelism (except specific tutorials on parallelism). However you can run most of the tutorial examples in parallel, see the topic on parallelism.
In case you work on your own PC or workstation, to make things easier, we suggest you define some handy environment variables by executing the following lines in the terminal:
export ABI_HOME=Replace_with_the_absolute_path_to_the_abinit_top_level_dir export PATH=$ABI_HOME/src/98_main/:$PATH export ABI_TESTS=$ABI_HOME/tests/ export ABI_PSPDIR=$ABI_TESTS/Psps_for_tests/ # Pseudopotentials used in examples.
Examples in this tutorial use these shell variables: copy and paste the code snippets into the terminal (remember to set ABI_HOME first!). The ‘export PATH’ line adds the directory containing the executables to your PATH so that you can invoke the code by simply typing abinit in the terminal instead of providing the absolute path.
1 General example of an almost converged GW calculation¶
Before beginning, you might consider to work in a different subdirectory as for the other tutorials. Why not Work_gw1?
At the end of tutorial 3, we computed the KS band structure of silicon. In this approximation, the band dispersion as well as the band widths are reasonable but the band gaps are qualitatively wrong. Now we will compute the band gaps much more accurately, using the so-called GW approximation.
We start by an example, in which we show how to perform in a single input file the calculation of the ground state density, the Kohn Sham band structure, the screening, and the GW corrections. We use reasonable values for the parameters of the calculation. The discussion on the convergence tests is postponed to the next paragraphs. We will see that GW calculations are much more time-consuming than the computation of the KS eigenvalues.
So, let us run immediately this calculation, and while it is running, we will explain what has been done.
cd $ABI_TESTS/tutorial/Input mkdir Work_gw1 cd Work_gw1 cp ../tgw1_x.files . # modify this file as usual (see tutorial 1) cp ../tgw1_1.in .
Then, issue:
abinit < tgw1_x.files > log 2> err &
Please run this job in background because it takes about 1 minute. In the meantime, you should read the following.
1.a The four steps of a GW calculation.¶
In order to perform a standard one-shot GW calculation one has to:
-
Run a converged Ground State calculation to obtain the self-consistent density.
-
Perform a non self-consistent run to compute the KS eigenvalues and the eigenfunctions including several empty states. Note that, unlike standard band structure calculations, here the KS states must be computed on a regular grid of k-points.
-
Use optdriver = 3 to compute the independent-particle susceptibility \chi^0 on a regular grid of q-points, for at least two frequencies (usually, \omega=0 and a large purely imaginary frequency - of the order of the plasmon frequency, a dozen of eV). The inverse dielectric matrix \epsilon^{-1} is then obtained via matrix inversion and stored in an external file (SCR). The list of q-points is automatically defined by the k-mesh used to generate the KS states in the previous step.
-
Use optdriver = 4 to compute the self-energy \Sigma matrix elements for a given set of k-points in order to obtain the GW quasiparticle energies. Note that the k-point must belong to the k-mesh used to generate the WFK file in step 2.
The flowchart diagram of a standard one-shot run is depicted in the figure below.
The input file tgw1_1.in has precisely that structure: there are four datasets.
The first dataset performs the SCF calculation to get the density. The second dataset reads the previous density file and performs a NSCF run including several empty states. The third dataset reads the WFK file produced in the previous step and drives the computation of susceptibility and dielectric matrices, producing another specialized file, tgw1_xo_DS2_SCR (_SCR for “Screening”, actually the inverse dielectric matrix \epsilon^{-1}). Then, in the fourth dataset, the code calculates the quasiparticle energies for the 4th and 5th bands at the \Gamma point.
So, you can edit this tgw1_1.in file.
../tgw1_x.in tgw1_x.out tgw1i tgw1o tgw1 ../../../Psps_for_tests/14si.pspnc
# Crystalline silicon # Calculation of the GW corrections # Dataset 1: ground state calculation to get the density # Dataset 2: NSCF run to produce the WFK file for 10 k-points in IBZ # Dataset 3: calculation of the screening (epsilon^-1 matrix for W) # Dataset 4: calculation of the Self-Energy matrix elements (GW corrections) ndtset 4 ############ # Dataset 1 ############ # SCF-GS run nband1 6 tolvrs1 1.0e-10 ############ # Dataset 2 ############ # Definition of parameters for the calculation of the WFK file nband2 40 # Number of (occ and empty) bands to be computed nbdbuf2 5 iscf2 -2 getden2 -1 tolwfr2 1.0d-18 # Will stop when this tolerance is achieved ############ # Dataset 3 ############ # Calculation of the screening (epsilon^-1 matrix) optdriver3 3 # Screening calculation getwfk3 -1 # Obtain WFK file from previous dataset nband3 17 # Bands to be used in the screening calculation ecuteps3 3.6 # Cut-off energy of the planewave set to represent the dielectric matrix. # It is important to adjust this parameter. ppmfrq3 16.7 eV # Imaginary frequency where to calculate the screening ############ # Dataset 4 ############ # Calculation of the Self-Energy matrix elements (GW corrections) optdriver4 4 # Self-Energy calculation getwfk4 -2 # Obtain WFK file from dataset 1 getscr4 -1 # Obtain SCR file from previous dataset nband4 30 # Bands to be used in the Self-Energy calculation ecutsigx4 8.0 # Dimension of the G sum in Sigma_x. # (the dimension in Sigma_c is controlled by ecuteps) nkptgw4 1 # number of k-point where to calculate the GW correction kptgw4 # k-points in reduced coordinates 0.000 0.000 0.000 bdgw4 4 5 # calculate GW corrections for bands from 4 to 5 # Data common to the three different datasets # Definition of the unit cell: fcc acell 3*10.217 # This is equivalent to 10.217 10.217 10.217 rprim 0.0 0.5 0.5 # FCC primitive vectors (to be scaled by acell) 0.5 0.0 0.5 0.5 0.5 0.0 # Definition of the atom types ntypat 1 # There is only one type of atom znucl 14 # The keyword "znucl" refers to the atomic number of the # possible type(s) of atom. The pseudopotential(s) # mentioned in the "files" file must correspond # to the type(s) of atom. Here, the only type is Silicon. # Definition of the atoms natom 2 # There are two atoms typat 1 1 # They both are of type 1, that is, Silicon. xred # Reduced coordinate of atoms 0.0 0.0 0.0 0.25 0.25 0.25 # Definition of the k-point grid ngkpt 2 2 2 nshiftk 4 shiftk 0.0 0.0 0.0 # These shifts will be the same for all grids 0.0 0.5 0.5 0.5 0.0 0.5 0.5 0.5 0.0 istwfk *1 # This is mandatory in all the GW steps. # Definition of the planewave basis set (at convergence 16 Rydberg 8 Hartree) ecut 8.0 # Maximal kinetic energy cut-off, in Hartree # Definition of the SCF procedure nstep 20 # Maximal number of SCF cycles #toldfe 1.0d-6 # Will stop when this tolerance is achieved on total energy diemac 12.0 # Although this is not mandatory, it is worth to # precondition the SCF cycle. The model dielectric # function used as the standard preconditioner # is described in the "dielng" input variable section. # Here, we follow the prescription for bulk silicon. #%%<BEGIN TEST_INFO> #%% [setup] #%% executable = abinit #%% [files] #%% files_to_test = #%% tgw1_1.out, tolnlines= 7, tolabs= 0.03, tolrel= 1.500e-01 #%% psp_files = 14si.pspnc #%% [paral_info] #%% max_nprocs = 4 #%% [extra_info] #%% authors = Unknown #%% keywords = GW #%% description = #%% Crystalline silicon #%% Calculation of the GW corrections #%% Dataset 1: ground state calculation and calculation of the WFK file for 10 k-points in IBZ #%% Dataset 2: calculation of the screening (epsilon^-1 matrix for W) #%% Dataset 3: calculation of the Self-Energy matrix elements (GW corrections) #%%<END TEST_INFO>
The dataset-independent part of this file (the last half of the file), contains the usual set of input variables describing the cell, atom types, number, position, planewave cut-off energy, SCF convergence parameters driving the KS band structure calculation. Then, for the fourth datasets, you will find specialized additional input variables.
1.b Generating the Kohn-Sham band structure: the WFK file.¶
Dataset 1 is a rather standard SCF calculation. It is worth noticing that we use tolvrs to stop the SCF cycle because we want a well-converged KS potential to be used in the subsequent NSCF calculation. Dataset 2 computes 40 bands and we set nbdbuf to 5 so that only the first 35 states must be converged within tolwfr. The 5 highest energy states are simply not considered when checking the convergence.
############ # Dataset 1 ############ # SCF-GS run nband1 6 tolvrs1 1.0e-10 ############ # Dataset 2 ############ # Definition of parameters for the calculation of the WFK file nband2 40 # Number of (occ and empty) bands to be computed nbdbuf2 5 iscf2 -2 getden2 -1 tolwfr2 1.0d-18 # Will stop when this tolerance is achieved
Important
The nbdbuf trick allows us to save several minimization steps because the last bands usually require more iterations to converge in the iterative diagonalization algorithms. Also note that it is a very good idea to increase significantly the value of nbdbuf when computing many empty states. As a rule of thumb use 10% of nband or even more in complicated systems. This can really make a huge difference at the level of the wall time.
1.c Generating the screening: the SCR file.¶
In dataset 3, the calculation of the screening (KS susceptibility \chi^0 and then inverse dielectric matrix \epsilon^{-1}) is performed. We need to set optdriver=3 to do that:
optdriver3 3 # Screening calculation
The getwfk input variable is similar to other “get” input variables of ABINIT:
getwfk3 -1 # Obtain WFK file from previous dataset
In this case, it tells the code to use the WFK file calculated in the previous dataset.
Then, three input variables describe the computation:
nband3 17 # Bands used in the screening calculation ecuteps3 3.6 # Cut-off energy of the planewave set to represent the dielectric matrix
In this case, we use 17 bands to calculate the KS response function \chi^{0}. The dimension of \chi^{0}, as well as all the other matrices (\chi, \epsilon^{-1}) is determined by the cut-off energy ecuteps = 3.6 Hartree, which yields 169 planewaves in our case.
Finally, we define the frequencies at which the screening must be evaluated: \omega=0.0 eV and the imaginary frequency \omega= i 16.7 eV. The latter is determined by the input variable ppmfrq
ppmfrq3 16.7 eV # Imaginary frequency where to calculate the screening
The two frequencies are used to calculate the plasmon-pole model parameters. For the non-zero frequency it is recommended to use a value close to the plasmon frequency for the plasmon-pole model to work well. Plasmons frequencies are usually close to 0.5 Hartree. The parameters for the screening calculation are not far from the ones that give converged Energy Loss Function (-\mathrm{Im} \epsilon^{-1}_{00}) spectra, so that one can start up by using indications from EELS calculations existing in literature.
1.d Computing the GW energies.¶
In dataset 4 the calculation of the Self-Energy matrix elements is performed. One needs to define the driver option as well as the _WFK and _SCR files.
optdriver4 4 # Self-Energy calculation getwfk4 -2 # Obtain WFK file from dataset 2 getscr4 -1 # Obtain SCR file from previous dataset
The getscr input variable is similar to other “get” input variables of ABINIT.
Then, comes the definition of parameters needed to compute the self-energy. As for the computation of the susceptibility and dielectric matrices, one must define the set of bands and two sets of planewaves:
nband4 30 # Bands to be used in the Self-Energy calculation ecutsigx4 8.0 # Dimension of the G sum in Sigma_x # (the dimension in Sigma_c is controlled by npweps)
In this case, nband controls the number of bands used to calculate the correlation part of the Self-Energy while ecutsigx gives the number of planewaves used to calculate \Sigma_x (the exchange part of the self-energy). The size of the planewave set used to compute \Sigma_c (the correlation part of the self-energy) is controlled by ecuteps and cannot be larger than the value used to generate the SCR file. For the initial convergence studies, it is advised to set ecutsigx to a value as high as ecut since, any way, this parameter is not much influential on the total computational time. Note that exact treatment of the exchange part requires, in principle, ecutsigx = 4 ecut.
Then, come the parameters defining the k-points and the band indices for which the quasiparticle energies will be computed:
nkptgw4 1 # number of k-point where to calculate the GW correction kptgw4 0.00 0.00 0.00 # k-points bdgw4 4 5 # calculate GW corrections for bands from 4 to 5
nkptgw defines the number of k-points for which the GW corrections will be computed. The k-point reduced coordinates are specified in kptgw. They must belong to the k-mesh used to generate the WFK file. Hence if you wish the GW correction in a particular k-point, you should choose a grid containing it. Usually this is done by taking the k-point grid where the convergence is achieved and shifting it such as at least one k-point is placed on the wished position in the Brillouin zone. bdgw gives the minimum/maximum band whose energies are calculated for each selected k-point.
There is an additional parameter, called zcut, related to the self-energy computation. It is meant to avoid some divergences that might occur in the calculation due to integrable poles along the integration path.
1.e Examination of the output file.¶
Let us hope that your calculation has been completed, and that we can examine the output file. Open tgw1_1.out in your preferred editor and find the section corresponding to DATASET 3.
.Version 9.0.0 of ABINIT .(MPI version, prepared for a x86_64_linux_gnu9.2 computer) .Copyright (C) 1998-2020 ABINIT group . ABINIT comes with ABSOLUTELY NO WARRANTY. It is free software, and you are welcome to redistribute it under certain conditions (GNU General Public License, see ~abinit/COPYING or http://www.gnu.org/copyleft/gpl.txt). ABINIT is a project of the Universite Catholique de Louvain, Corning Inc. and other collaborators, see ~abinit/doc/developers/contributors.txt . Please read https://docs.abinit.org/theory/acknowledgments for suggested acknowledgments of the ABINIT effort. For more information, see https://www.abinit.org . .Starting date : Mon 24 Feb 2020. - ( at 16h45 ) - input file -> /home/gmatteo/git_repos/abinit/_abiref_gnu9.2_openmpi/tests/Test_suite/tutorial_tgw1_1/tgw1_1.in - output file -> tgw1_1.out - root for input files -> tgw1_1i - root for output files -> tgw1_1o DATASET 1 : space group Fd -3 m (#227); Bravais cF (face-center cubic) ================================================================================ Values of the parameters that define the memory need for DATASET 1. intxc = 0 ionmov = 0 iscf = 7 lmnmax = 2 lnmax = 2 mgfft = 20 mpssoang = 3 mqgrid = 3001 natom = 2 nloc_mem = 1 nspden = 1 nspinor = 1 nsppol = 1 nsym = 48 n1xccc = 2501 ntypat = 1 occopt = 1 xclevel = 1 - mband = 6 mffmem = 1 mkmem = 6 mpw = 303 nfft = 8000 nkpt = 6 ================================================================================ P This job should need less than 3.014 Mbytes of memory. Rough estimation (10% accuracy) of disk space for files : _ WF disk file : 0.168 Mbytes ; DEN or POT disk file : 0.063 Mbytes. ================================================================================ DATASET 2 : space group Fd -3 m (#227); Bravais cF (face-center cubic) ================================================================================ Values of the parameters that define the memory need for DATASET 2. intxc = 0 ionmov = 0 iscf = -2 lmnmax = 2 lnmax = 2 mgfft = 20 mpssoang = 3 mqgrid = 3001 natom = 2 nloc_mem = 1 nspden = 1 nspinor = 1 nsppol = 1 nsym = 48 n1xccc = 2501 ntypat = 1 occopt = 1 xclevel = 1 - mband = 40 mffmem = 1 mkmem = 6 mpw = 303 nfft = 8000 nkpt = 6 ================================================================================ P This job should need less than 3.053 Mbytes of memory. Rough estimation (10% accuracy) of disk space for files : _ WF disk file : 1.112 Mbytes ; DEN or POT disk file : 0.063 Mbytes. ================================================================================ DATASET 3 : space group Fd -3 m (#227); Bravais cF (face-center cubic) ================================================================================ Values of the parameters that define the memory need for DATASET 3. intxc = 0 ionmov = 0 iscf = 7 lmnmax = 2 lnmax = 2 mgfft = 18 mpssoang = 3 mqgrid = 3001 natom = 2 nloc_mem = 1 nspden = 1 nspinor = 1 nsppol = 1 nsym = 48 n1xccc = 2501 ntypat = 1 occopt = 1 xclevel = 1 - mband = 17 mffmem = 1 mkmem = 6 mpw = 266 nfft = 5832 nkpt = 6 ================================================================================ P This job should need less than 2.669 Mbytes of memory. Rough estimation (10% accuracy) of disk space for files : _ WF disk file : 0.416 Mbytes ; DEN or POT disk file : 0.046 Mbytes. ================================================================================ DATASET 4 : space group Fd -3 m (#227); Bravais cF (face-center cubic) ================================================================================ Values of the parameters that define the memory need for DATASET 4. intxc = 0 ionmov = 0 iscf = 7 lmnmax = 2 lnmax = 2 mgfft = 18 mpssoang = 3 mqgrid = 3001 natom = 2 nloc_mem = 1 nspden = 1 nspinor = 1 nsppol = 1 nsym = 48 n1xccc = 2501 ntypat = 1 occopt = 1 xclevel = 1 - mband = 30 mffmem = 1 mkmem = 6 mpw = 266 nfft = 5832 nkpt = 6 ================================================================================ P This job should need less than 3.013 Mbytes of memory. Rough estimation (10% accuracy) of disk space for files : _ WF disk file : 0.733 Mbytes ; DEN or POT disk file : 0.046 Mbytes. ================================================================================ -------------------------------------------------------------------------------- ------------- Echo of variables that govern the present computation ------------ -------------------------------------------------------------------------------- - - outvars: echo of selected default values - iomode0 = 0 , fftalg0 =312 , wfoptalg0 = 0 - - outvars: echo of global parameters not present in the input file - max_nthreads = 0 - -outvars: echo values of preprocessed input variables -------- acell 1.0217000000E+01 1.0217000000E+01 1.0217000000E+01 Bohr amu 2.80855000E+01 bdgw4 4 5 diemac 1.20000000E+01 ecut1 8.00000000E+00 Hartree ecut2 8.00000000E+00 Hartree ecut3 7.56385066E+00 Hartree ecut4 7.56385066E+00 Hartree ecuteps1 0.00000000E+00 Hartree ecuteps2 0.00000000E+00 Hartree ecuteps3 3.59282906E+00 Hartree ecuteps4 0.00000000E+00 Hartree ecutsigx1 0.00000000E+00 Hartree ecutsigx2 0.00000000E+00 Hartree ecutsigx3 0.00000000E+00 Hartree ecutsigx4 7.56385066E+00 Hartree ecutwfn1 0.00000000E+00 Hartree ecutwfn2 0.00000000E+00 Hartree ecutwfn3 7.56385066E+00 Hartree ecutwfn4 7.56385066E+00 Hartree - fftalg 312 getden1 0 getden2 -1 getden3 0 getden4 0 getscr1 0 getscr2 0 getscr3 0 getscr4 -1 getwfk1 0 getwfk2 0 getwfk3 -1 getwfk4 -2 iscf1 7 iscf2 -2 iscf3 7 iscf4 7 istwfk 0 0 1 0 1 1 jdtset 1 2 3 4 kpt -2.50000000E-01 -2.50000000E-01 0.00000000E+00 -2.50000000E-01 2.50000000E-01 0.00000000E+00 5.00000000E-01 5.00000000E-01 0.00000000E+00 -2.50000000E-01 5.00000000E-01 2.50000000E-01 5.00000000E-01 0.00000000E+00 0.00000000E+00 0.00000000E+00 0.00000000E+00 0.00000000E+00 kptgw4 0.00000000E+00 0.00000000E+00 0.00000000E+00 kptrlatt 2 -2 2 -2 2 2 -2 -2 2 kptrlen 2.04340000E+01 P mkmem 6 natom 2 nband1 6 nband2 40 nband3 17 nband4 30 nbdbuf1 0 nbdbuf2 5 nbdbuf3 0 nbdbuf4 0 ndtset 4 ngfft1 20 20 20 ngfft2 20 20 20 ngfft3 18 18 18 ngfft4 18 18 18 nkpt 6 nkptgw1 0 nkptgw2 0 nkptgw3 0 nkptgw4 1 npweps1 0 npweps2 0 npweps3 89 npweps4 0 npwsigx1 0 npwsigx2 0 npwsigx3 0 npwsigx4 283 npwwfn1 0 npwwfn2 0 npwwfn3 283 npwwfn4 283 nstep 20 nsym 48 ntypat 1 occ1 2.000000 2.000000 2.000000 2.000000 0.000000 0.000000 occ3 2.000000 2.000000 2.000000 2.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 occ4 2.000000 2.000000 2.000000 2.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 optdriver1 0 optdriver2 0 optdriver3 3 optdriver4 4 ppmfrq1 0.00000000E+00 Hartree ppmfrq2 0.00000000E+00 Hartree ppmfrq3 6.13713734E-01 Hartree ppmfrq4 0.00000000E+00 Hartree rprim 0.0000000000E+00 5.0000000000E-01 5.0000000000E-01 5.0000000000E-01 0.0000000000E+00 5.0000000000E-01 5.0000000000E-01 5.0000000000E-01 0.0000000000E+00 spgroup 227 symrel 1 0 0 0 1 0 0 0 1 -1 0 0 0 -1 0 0 0 -1 0 -1 1 0 -1 0 1 -1 0 0 1 -1 0 1 0 -1 1 0 -1 0 0 -1 0 1 -1 1 0 1 0 0 1 0 -1 1 -1 0 0 1 -1 1 0 -1 0 0 -1 0 -1 1 -1 0 1 0 0 1 -1 0 0 -1 1 0 -1 0 1 1 0 0 1 -1 0 1 0 -1 0 -1 1 1 -1 0 0 -1 0 0 1 -1 -1 1 0 0 1 0 1 0 0 0 0 1 0 1 0 -1 0 0 0 0 -1 0 -1 0 0 1 -1 0 0 -1 1 0 -1 0 -1 1 0 0 1 -1 0 1 -1 0 1 -1 1 0 -1 0 0 1 0 -1 1 -1 0 1 0 0 0 -1 0 1 -1 0 0 -1 1 0 1 0 -1 1 0 0 1 -1 1 0 -1 0 0 -1 0 1 -1 -1 0 1 0 0 1 0 -1 1 0 1 0 0 0 1 1 0 0 0 -1 0 0 0 -1 -1 0 0 1 0 -1 0 1 -1 0 0 -1 -1 0 1 0 -1 1 0 0 1 0 -1 0 0 -1 1 1 -1 0 0 1 0 0 1 -1 -1 1 0 -1 0 1 -1 0 0 -1 1 0 1 0 -1 1 0 0 1 -1 0 0 1 0 1 0 0 0 0 1 0 -1 0 -1 0 0 0 0 -1 0 0 -1 0 1 -1 1 0 -1 0 0 1 0 -1 1 -1 0 1 1 -1 0 0 -1 1 0 -1 0 -1 1 0 0 1 -1 0 1 0 0 0 1 1 0 0 0 1 0 0 0 -1 -1 0 0 0 -1 0 -1 1 0 -1 0 0 -1 0 1 1 -1 0 1 0 0 1 0 -1 0 0 1 0 1 0 1 0 0 0 0 -1 0 -1 0 -1 0 0 1 -1 0 0 -1 0 0 -1 1 -1 1 0 0 1 0 0 1 -1 0 0 -1 1 0 -1 0 1 -1 0 0 1 -1 0 1 0 -1 1 -1 1 0 -1 0 1 -1 0 0 1 -1 0 1 0 -1 1 0 0 tnons 0.0000000 0.0000000 0.0000000 0.2500000 0.2500000 0.2500000 0.0000000 0.0000000 0.0000000 0.2500000 0.2500000 0.2500000 0.0000000 0.0000000 0.0000000 0.2500000 0.2500000 0.2500000 0.0000000 0.0000000 0.0000000 0.2500000 0.2500000 0.2500000 0.0000000 0.0000000 0.0000000 0.2500000 0.2500000 0.2500000 0.0000000 0.0000000 0.0000000 0.2500000 0.2500000 0.2500000 0.0000000 0.0000000 0.0000000 0.2500000 0.2500000 0.2500000 0.0000000 0.0000000 0.0000000 0.2500000 0.2500000 0.2500000 0.0000000 0.0000000 0.0000000 0.2500000 0.2500000 0.2500000 0.0000000 0.0000000 0.0000000 0.2500000 0.2500000 0.2500000 0.0000000 0.0000000 0.0000000 0.2500000 0.2500000 0.2500000 0.0000000 0.0000000 0.0000000 0.2500000 0.2500000 0.2500000 0.0000000 0.0000000 0.0000000 0.2500000 0.2500000 0.2500000 0.0000000 0.0000000 0.0000000 0.2500000 0.2500000 0.2500000 0.0000000 0.0000000 0.0000000 0.2500000 0.2500000 0.2500000 0.0000000 0.0000000 0.0000000 0.2500000 0.2500000 0.2500000 0.0000000 0.0000000 0.0000000 0.2500000 0.2500000 0.2500000 0.0000000 0.0000000 0.0000000 0.2500000 0.2500000 0.2500000 0.0000000 0.0000000 0.0000000 0.2500000 0.2500000 0.2500000 0.0000000 0.0000000 0.0000000 0.2500000 0.2500000 0.2500000 0.0000000 0.0000000 0.0000000 0.2500000 0.2500000 0.2500000 0.0000000 0.0000000 0.0000000 0.2500000 0.2500000 0.2500000 0.0000000 0.0000000 0.0000000 0.2500000 0.2500000 0.2500000 0.0000000 0.0000000 0.0000000 0.2500000 0.2500000 0.2500000 tolvrs1 1.00000000E-10 tolvrs2 0.00000000E+00 tolvrs3 0.00000000E+00 tolvrs4 0.00000000E+00 tolwfr1 0.00000000E+00 tolwfr2 1.00000000E-18 tolwfr3 0.00000000E+00 tolwfr4 0.00000000E+00 typat 1 1 wtk 0.18750 0.37500 0.09375 0.18750 0.12500 0.03125 xangst 0.0000000000E+00 0.0000000000E+00 0.0000000000E+00 1.3516508850E+00 1.3516508850E+00 1.3516508850E+00 xcart 0.0000000000E+00 0.0000000000E+00 0.0000000000E+00 2.5542500000E+00 2.5542500000E+00 2.5542500000E+00 xred 0.0000000000E+00 0.0000000000E+00 0.0000000000E+00 2.5000000000E-01 2.5000000000E-01 2.5000000000E-01 znucl 14.00000 ================================================================================ chkinp: Checking input parameters for consistency, jdtset= 1. chkinp: Checking input parameters for consistency, jdtset= 2. chkinp: Checking input parameters for consistency, jdtset= 3. chkinp: Checking input parameters for consistency, jdtset= 4. ================================================================================ == DATASET 1 ================================================================== - mpi_nproc: 1, omp_nthreads: -1 (-1 if OMP is not activated) --- !DatasetInfo iteration_state: {dtset: 1, } dimensions: {natom: 2, nkpt: 6, mband: 6, nsppol: 1, nspinor: 1, nspden: 1, mpw: 303, } cutoff_energies: {ecut: 8.0, pawecutdg: -1.0, } electrons: {nelect: 8.00000000E+00, charge: 0.00000000E+00, occopt: 1.00000000E+00, tsmear: 1.00000000E-02, } meta: {optdriver: 0, ionmov: 0, optcell: 0, iscf: 7, paral_kgb: 0, } ... Exchange-correlation functional for the present dataset will be: LDA: new Teter (4/93) with spin-polarized option - ixc=1 Citation for XC functional: S. Goedecker, M. Teter, J. Huetter, PRB 54, 1703 (1996) Real(R)+Recip(G) space primitive vectors, cartesian coordinates (Bohr,Bohr^-1): R(1)= 0.0000000 5.1085000 5.1085000 G(1)= -0.0978761 0.0978761 0.0978761 R(2)= 5.1085000 0.0000000 5.1085000 G(2)= 0.0978761 -0.0978761 0.0978761 R(3)= 5.1085000 5.1085000 0.0000000 G(3)= 0.0978761 0.0978761 -0.0978761 Unit cell volume ucvol= 2.6663072E+02 bohr^3 Angles (23,13,12)= 6.00000000E+01 6.00000000E+01 6.00000000E+01 degrees getcut: wavevector= 0.0000 0.0000 0.0000 ngfft= 20 20 20 ecut(hartree)= 8.000 => boxcut(ratio)= 2.17426 --- Pseudopotential description ------------------------------------------------ - pspini: atom type 1 psp file is /home/gmatteo/git_repos/abinit/tests/Psps_for_tests/14si.pspnc - pspatm: opening atomic psp file /home/gmatteo/git_repos/abinit/tests/Psps_for_tests/14si.pspnc - Troullier-Martins psp for element Si Thu Oct 27 17:31:21 EDT 1994 - 14.00000 4.00000 940714 znucl, zion, pspdat 1 1 2 2 2001 0.00000 pspcod,pspxc,lmax,lloc,mmax,r2well 0 5.907 14.692 1 2.0872718 l,e99.0,e99.9,nproj,rcpsp 0.00000000 0.00000000 0.00000000 0.00000000 rms, ekb1, ekb2, epsatm 1 2.617 4.181 1 2.0872718 l,e99.0,e99.9,nproj,rcpsp 0.00000000 0.00000000 0.00000000 0.00000000 rms, ekb1, ekb2, epsatm 2 0.000 0.000 0 2.0872718 l,e99.0,e99.9,nproj,rcpsp 0.00000000 0.00000000 0.00000000 0.00000000 rms, ekb1, ekb2, epsatm 1.80626423934776 0.22824404341771 1.17378968127746 rchrg,fchrg,qchrg pspatm : epsatm= 1.43386982 --- l ekb(1:nproj) --> 0 3.287949 1 1.849886 pspatm: atomic psp has been read and splines computed 2.29419171E+01 ecore*ucvol(ha*bohr**3) -------------------------------------------------------------------------------- _setup2: Arith. and geom. avg. npw (full set) are 286.594 286.400 ================================================================================ --- !BeginCycle iteration_state: {dtset: 1, } solver: {iscf: 7, nstep: 20, nline: 4, wfoptalg: 0, } tolerances: {tolvrs: 1.00E-10, } ... iter Etot(hartree) deltaE(h) residm vres2 ETOT 1 -8.8582779584510 -8.858E+00 6.600E-03 6.417E+00 ETOT 2 -8.8641756264350 -5.898E-03 4.786E-05 1.710E-01 ETOT 3 -8.8642530622065 -7.744E-05 7.462E-06 3.886E-03 ETOT 4 -8.8642542690300 -1.207E-06 3.176E-06 1.275E-05 ETOT 5 -8.8642542741682 -5.138E-09 4.336E-07 6.761E-09 ETOT 6 -8.8642542741704 -2.252E-12 2.993E-07 3.096E-11 At SCF step 6 vres2 = 3.10E-11 < tolvrs= 1.00E-10 =>converged. Cartesian components of stress tensor (hartree/bohr^3) sigma(1 1)= 9.87319501E-07 sigma(3 2)= 0.00000000E+00 sigma(2 2)= 9.87319501E-07 sigma(3 1)= 0.00000000E+00 sigma(3 3)= 9.87319501E-07 sigma(2 1)= 0.00000000E+00 --- !ResultsGS iteration_state: {dtset: 1, } comment : Summary of ground state results lattice_vectors: - [ 0.0000000, 5.1085000, 5.1085000, ] - [ 5.1085000, 0.0000000, 5.1085000, ] - [ 5.1085000, 5.1085000, 0.0000000, ] lattice_lengths: [ 7.22451, 7.22451, 7.22451, ] lattice_angles: [ 60.000, 60.000, 60.000, ] # degrees, (23, 13, 12) lattice_volume: 2.6663072E+02 convergence: {deltae: -2.252E-12, res2: 3.096E-11, residm: 2.993E-07, diffor: null, } etotal : -8.86425427E+00 entropy : 0.00000000E+00 fermie : 2.19278031E-01 cartesian_stress_tensor: # hartree/bohr^3 - [ 9.87319501E-07, 0.00000000E+00, 0.00000000E+00, ] - [ 0.00000000E+00, 9.87319501E-07, 0.00000000E+00, ] - [ 0.00000000E+00, 0.00000000E+00, 9.87319501E-07, ] pressure_GPa: -2.9048E-02 xred : - [ 0.0000E+00, 0.0000E+00, 0.0000E+00, Si] - [ 2.5000E-01, 2.5000E-01, 2.5000E-01, Si] cartesian_forces: # hartree/bohr - [ -0.00000000E+00, -0.00000000E+00, -0.00000000E+00, ] - [ -0.00000000E+00, -0.00000000E+00, -0.00000000E+00, ] force_length_stats: {min: 0.00000000E+00, max: 0.00000000E+00, mean: 0.00000000E+00, } ... Integrated electronic density in atomic spheres: ------------------------------------------------ Atom Sphere_radius Integrated_density 1 2.00000 1.72095363 2 2.00000 1.72095363 ================================================================================ ----iterations are completed or convergence reached---- Mean square residual over all n,k,spin= 84.662E-10; max= 29.932E-08 reduced coordinates (array xred) for 2 atoms 0.000000000000 0.000000000000 0.000000000000 0.250000000000 0.250000000000 0.250000000000 rms dE/dt= 0.0000E+00; max dE/dt= 0.0000E+00; dE/dt below (all hartree) 1 0.000000000000 0.000000000000 0.000000000000 2 0.000000000000 0.000000000000 0.000000000000 cartesian coordinates (angstrom) at end: 1 0.00000000000000 0.00000000000000 0.00000000000000 2 1.35165088504101 1.35165088504101 1.35165088504101 cartesian forces (hartree/bohr) at end: 1 -0.00000000000000 -0.00000000000000 -0.00000000000000 2 -0.00000000000000 -0.00000000000000 -0.00000000000000 frms,max,avg= 0.0000000E+00 0.0000000E+00 0.000E+00 0.000E+00 0.000E+00 h/b cartesian forces (eV/Angstrom) at end: 1 -0.00000000000000 -0.00000000000000 -0.00000000000000 2 -0.00000000000000 -0.00000000000000 -0.00000000000000 frms,max,avg= 0.0000000E+00 0.0000000E+00 0.000E+00 0.000E+00 0.000E+00 e/A length scales= 10.217000000000 10.217000000000 10.217000000000 bohr = 5.406603540164 5.406603540164 5.406603540164 angstroms prteigrs : about to open file tgw1_1o_DS1_EIG Fermi (or HOMO) energy (hartree) = 0.21928 Average Vxc (hartree)= -0.35388 Eigenvalues (hartree) for nkpt= 6 k points: kpt# 1, nband= 6, wtk= 0.18750, kpt= -0.2500 -0.2500 0.0000 (reduced coord) -0.18502 0.08915 0.14773 0.14773 0.25608 0.33012 prteigrs : prtvol=0 or 1, do not print more k-points. --- !EnergyTerms iteration_state : {dtset: 1, } comment : Components of total free energy in Hartree kinetic : 3.05324036764975E+00 hartree : 5.51420121188750E-01 xc : -3.54423040234180E+00 Ewald energy : -8.43581958561899E+00 psp_core : 8.60437873155177E-02 local_psp : -2.45083905740684E+00 non_local_psp : 1.87593049504319E+00 total_energy : -8.86425427417042E+00 total_energy_eV : -2.41208625687097E+02 band_energy : 2.51249507604680E-01 ... Cartesian components of stress tensor (hartree/bohr^3) sigma(1 1)= 9.87319501E-07 sigma(3 2)= 0.00000000E+00 sigma(2 2)= 9.87319501E-07 sigma(3 1)= 0.00000000E+00 sigma(3 3)= 9.87319501E-07 sigma(2 1)= 0.00000000E+00 -Cartesian components of stress tensor (GPa) [Pressure= -2.9048E-02 GPa] - sigma(1 1)= 2.90479377E-02 sigma(3 2)= 0.00000000E+00 - sigma(2 2)= 2.90479377E-02 sigma(3 1)= 0.00000000E+00 - sigma(3 3)= 2.90479377E-02 sigma(2 1)= 0.00000000E+00 ================================================================================ == DATASET 2 ================================================================== - mpi_nproc: 1, omp_nthreads: -1 (-1 if OMP is not activated) --- !DatasetInfo iteration_state: {dtset: 2, } dimensions: {natom: 2, nkpt: 6, mband: 40, nsppol: 1, nspinor: 1, nspden: 1, mpw: 303, } cutoff_energies: {ecut: 8.0, pawecutdg: -1.0, } electrons: {nelect: 8.00000000E+00, charge: 0.00000000E+00, occopt: 1.00000000E+00, tsmear: 1.00000000E-02, } meta: {optdriver: 0, ionmov: 0, optcell: 0, iscf: -2, paral_kgb: 0, } ... mkfilename : getden/=0, take file _DEN from output of DATASET 1. Exchange-correlation functional for the present dataset will be: LDA: new Teter (4/93) with spin-polarized option - ixc=1 Citation for XC functional: S. Goedecker, M. Teter, J. Huetter, PRB 54, 1703 (1996) Real(R)+Recip(G) space primitive vectors, cartesian coordinates (Bohr,Bohr^-1): R(1)= 0.0000000 5.1085000 5.1085000 G(1)= -0.0978761 0.0978761 0.0978761 R(2)= 5.1085000 0.0000000 5.1085000 G(2)= 0.0978761 -0.0978761 0.0978761 R(3)= 5.1085000 5.1085000 0.0000000 G(3)= 0.0978761 0.0978761 -0.0978761 Unit cell volume ucvol= 2.6663072E+02 bohr^3 Angles (23,13,12)= 6.00000000E+01 6.00000000E+01 6.00000000E+01 degrees getcut: wavevector= 0.0000 0.0000 0.0000 ngfft= 20 20 20 ecut(hartree)= 8.000 => boxcut(ratio)= 2.17426 -------------------------------------------------------------------------------- ================================================================================ prteigrs : about to open file tgw1_1o_DS2_EIG Non-SCF case, kpt 1 ( -0.25000 -0.25000 0.00000), residuals and eigenvalues= 1.82E-19 1.09E-19 6.71E-20 1.63E-19 7.84E-20 8.47E-19 5.08E-20 5.86E-20 3.79E-20 6.55E-20 2.17E-19 3.62E-20 7.98E-19 3.19E-19 2.32E-19 1.26E-19 3.63E-19 7.16E-19 1.16E-19 9.23E-19 8.24E-20 2.78E-19 4.69E-19 7.39E-19 6.01E-19 8.92E-20 5.85E-20 8.96E-19 7.87E-19 3.08E-19 3.24E-19 4.88E-19 6.30E-19 7.75E-19 8.05E-19 3.02E-17 8.27E-14 8.56E-14 7.04E-11 7.32E-07 -1.8502E-01 8.9150E-02 1.4773E-01 1.4773E-01 2.5608E-01 3.3012E-01 4.3046E-01 4.3046E-01 5.1148E-01 5.5223E-01 6.2983E-01 6.9695E-01 6.9695E-01 7.2030E-01 8.5843E-01 8.5843E-01 9.0176E-01 9.5958E-01 9.8832E-01 1.1388E+00 1.1880E+00 1.2244E+00 1.2244E+00 1.3032E+00 1.3032E+00 1.3173E+00 1.3273E+00 1.4938E+00 1.5078E+00 1.5334E+00 1.5334E+00 1.5944E+00 1.5944E+00 1.5966E+00 1.6450E+00 1.6671E+00 1.7256E+00 1.7256E+00 1.7528E+00 1.8311E+00 prteigrs : nnsclo,ikpt= 20 1 max resid (incl. the buffer)= 7.32068E-07 prteigrs : prtvol=0 or 1, do not print more k-points. prteigrs : nnsclo,ikpt= 20 2 max resid (incl. the buffer)= 3.39096E-06 prteigrs : nnsclo,ikpt= 20 3 max resid (incl. the buffer)= 3.54667E-16 prteigrs : nnsclo,ikpt= 20 4 max resid (incl. the buffer)= 1.15277E-10 prteigrs : nnsclo,ikpt= 20 5 max resid (incl. the buffer)= 1.30604E-13 prteigrs : nnsclo,ikpt= 20 6 max resid (incl. the buffer)= 4.13457E-06 --- !ResultsGS iteration_state: {dtset: 2, } comment : Summary of ground state results lattice_vectors: - [ 0.0000000, 5.1085000, 5.1085000, ] - [ 5.1085000, 0.0000000, 5.1085000, ] - [ 5.1085000, 5.1085000, 0.0000000, ] lattice_lengths: [ 7.22451, 7.22451, 7.22451, ] lattice_angles: [ 60.000, 60.000, 60.000, ] # degrees, (23, 13, 12) lattice_volume: 2.6663072E+02 convergence: {deltae: 0.000E+00, res2: 0.000E+00, residm: 9.789E-19, diffor: 0.000E+00, } etotal : -8.86425427E+00 entropy : 0.00000000E+00 fermie : 2.19278031E-01 cartesian_stress_tensor: null pressure_GPa: null xred : - [ 0.0000E+00, 0.0000E+00, 0.0000E+00, Si] - [ 2.5000E-01, 2.5000E-01, 2.5000E-01, Si] cartesian_forces: null force_length_stats: {min: null, max: null, mean: null, } ... Integrated electronic density in atomic spheres: ------------------------------------------------ Atom Sphere_radius Integrated_density 1 2.00000 1.72095363 2 2.00000 1.72095363 ================================================================================ ----iterations are completed or convergence reached---- Mean square residual over all n,k,spin= 33.755E-20; max= 97.890E-20 reduced coordinates (array xred) for 2 atoms 0.000000000000 0.000000000000 0.000000000000 0.250000000000 0.250000000000 0.250000000000 cartesian coordinates (angstrom) at end: 1 0.00000000000000 0.00000000000000 0.00000000000000 2 1.35165088504101 1.35165088504101 1.35165088504101 length scales= 10.217000000000 10.217000000000 10.217000000000 bohr = 5.406603540164 5.406603540164 5.406603540164 angstroms prteigrs : about to open file tgw1_1o_DS2_EIG Eigenvalues (hartree) for nkpt= 6 k points: kpt# 1, nband= 40, wtk= 0.18750, kpt= -0.2500 -0.2500 0.0000 (reduced coord) -0.18502 0.08915 0.14773 0.14773 0.25608 0.33012 0.43046 0.43046 0.51148 0.55223 0.62983 0.69695 0.69695 0.72030 0.85843 0.85843 0.90176 0.95958 0.98832 1.13879 1.18801 1.22440 1.22440 1.30322 1.30322 1.31734 1.32731 1.49384 1.50783 1.53336 1.53336 1.59442 1.59442 1.59659 1.64504 1.66705 1.72558 1.72558 1.75279 1.83111 prteigrs : prtvol=0 or 1, do not print more k-points. ================================================================================ == DATASET 3 ================================================================== - mpi_nproc: 1, omp_nthreads: -1 (-1 if OMP is not activated) --- !DatasetInfo iteration_state: {dtset: 3, } dimensions: {natom: 2, nkpt: 6, mband: 17, nsppol: 1, nspinor: 1, nspden: 1, mpw: 266, } cutoff_energies: {ecut: 7.6, pawecutdg: -1.0, } electrons: {nelect: 8.00000000E+00, charge: 0.00000000E+00, occopt: 1.00000000E+00, tsmear: 1.00000000E-02, } meta: {optdriver: 3, gwcalctyp: 0, } ... mkfilename : getwfk/=0, take file _WFK from output of DATASET 2. Exchange-correlation functional for the present dataset will be: LDA: new Teter (4/93) with spin-polarized option - ixc=1 Citation for XC functional: S. Goedecker, M. Teter, J. Huetter, PRB 54, 1703 (1996) SCREENING: Calculation of the susceptibility and dielectric matrices Based on a program developped by R.W. Godby, V. Olevano, G. Onida, and L. Reining. Incorporated in ABINIT by V. Olevano, G.-M. Rignanese, and M. Torrent. .Using double precision arithmetic ; gwpc = 8 Real(R)+Recip(G) space primitive vectors, cartesian coordinates (Bohr,Bohr^-1): R(1)= 0.0000000 5.1085000 5.1085000 G(1)= -0.0978761 0.0978761 0.0978761 R(2)= 5.1085000 0.0000000 5.1085000 G(2)= 0.0978761 -0.0978761 0.0978761 R(3)= 5.1085000 5.1085000 0.0000000 G(3)= 0.0978761 0.0978761 -0.0978761 Unit cell volume ucvol= 2.6663072E+02 bohr^3 Angles (23,13,12)= 6.00000000E+01 6.00000000E+01 6.00000000E+01 degrees --- Pseudopotential description ------------------------------------------------ - pspini: atom type 1 psp file is /home/gmatteo/git_repos/abinit/tests/Psps_for_tests/14si.pspnc - pspatm: opening atomic psp file /home/gmatteo/git_repos/abinit/tests/Psps_for_tests/14si.pspnc - Troullier-Martins psp for element Si Thu Oct 27 17:31:21 EDT 1994 - 14.00000 4.00000 940714 znucl, zion, pspdat 1 1 2 2 2001 0.00000 pspcod,pspxc,lmax,lloc,mmax,r2well 0 5.907 14.692 1 2.0872718 l,e99.0,e99.9,nproj,rcpsp 0.00000000 0.00000000 0.00000000 0.00000000 rms, ekb1, ekb2, epsatm 1 2.617 4.181 1 2.0872718 l,e99.0,e99.9,nproj,rcpsp 0.00000000 0.00000000 0.00000000 0.00000000 rms, ekb1, ekb2, epsatm 2 0.000 0.000 0 2.0872718 l,e99.0,e99.9,nproj,rcpsp 0.00000000 0.00000000 0.00000000 0.00000000 rms, ekb1, ekb2, epsatm 1.80626423934776 0.22824404341771 1.17378968127746 rchrg,fchrg,qchrg pspatm : epsatm= 1.43386982 --- l ekb(1:nproj) --> 0 3.287949 1 1.849886 pspatm: atomic psp has been read and splines computed -------------------------------------------------------------------------------- ==== K-mesh for the wavefunctions ==== Number of points in the irreducible wedge : 6 Reduced coordinates and weights : 1) -2.50000000E-01 -2.50000000E-01 0.00000000E+00 0.18750 2) -2.50000000E-01 2.50000000E-01 0.00000000E+00 0.37500 3) 5.00000000E-01 5.00000000E-01 0.00000000E+00 0.09375 4) -2.50000000E-01 5.00000000E-01 2.50000000E-01 0.18750 5) 5.00000000E-01 0.00000000E+00 0.00000000E+00 0.12500 6) 0.00000000E+00 0.00000000E+00 0.00000000E+00 0.03125 Together with 48 symmetry operations and time-reversal symmetry yields 32 points in the full Brillouin Zone. ==== Q-mesh for the screening function ==== Number of points in the irreducible wedge : 6 Reduced coordinates and weights : 1) 0.00000000E+00 0.00000000E+00 0.00000000E+00 0.03125 2) 5.00000000E-01 5.00000000E-01 0.00000000E+00 0.09375 3) 5.00000000E-01 2.50000000E-01 2.50000000E-01 0.37500 4) 0.00000000E+00 5.00000000E-01 0.00000000E+00 0.12500 5) 5.00000000E-01 -2.50000000E-01 2.50000000E-01 0.18750 6) 0.00000000E+00 -2.50000000E-01 -2.50000000E-01 0.18750 Together with 48 symmetry operations and time-reversal symmetry yields 32 points in the full Brillouin Zone. setmesh: FFT mesh size selected = 20x 20x 20 total number of points = 8000 - screening: taking advantage of time-reversal symmetry - Maximum band index for partially occupied states nbvw = 4 - Remaining bands to be divided among processors nbcw = 13 - Number of bands treated by each node ~13 Number of electrons calculated from density = 7.9999; Expected = 8.0000 average of density, n = 0.030004 r_s = 1.9964 omega_plasma = 16.7087 [eV] calculating chi0 at frequencies [eV] : 1 0.000000E+00 0.000000E+00 2 0.000000E+00 1.670000E+01 -------------------------------------------------------------------------------- q-point number 1 q = ( 0.000000, 0.000000, 0.000000) [r.l.u.] -------------------------------------------------------------------------------- chi0(G,G') at the 1 th omega 0.0000 0.0000 [eV] 1 2 3 4 5 6 7 8 9 1 -0.000 0.000 -0.000 -0.000 0.000 0.000 -0.000 -0.000 0.000 -0.000 0.000 0.000 -0.000 -0.000 0.000 0.000 -0.000 -0.000 2 0.000 -16.032 -0.000 -0.025 0.000 -0.021 0.000 -0.025 -0.000 -0.000 -0.000 -5.108 0.000 -0.282 -0.000 -0.252 0.000 -0.256 chi0(G,G') at the 2 th omega 0.0000 16.7000 [eV] 1 2 3 4 5 6 7 8 9 1 -0.000 0.000 -0.000 -0.000 0.000 0.000 -0.000 -0.000 0.000 -0.000 0.000 0.000 -0.000 -0.000 0.000 0.000 -0.000 -0.000 2 0.000 -5.610 -0.000 0.068 0.000 0.070 0.000 0.067 -0.000 -0.000 -0.000 -1.068 0.000 -0.119 -0.000 -0.099 0.000 -0.102 For q-point: 0.000010 0.000020 0.000030 dielectric constant = 24.4980 dielectric constant without local fields = 27.1361 Average fulfillment of the sum rule on Im[epsilon] for q-point 1 : 37.14 [%] Heads and wings of the symmetrical epsilon^-1(G,G') Upper and lower wings at the 1 th omega 0.0000 0.0000 [eV] 1 2 3 4 5 6 7 8 9 0.041 0.004 -0.004 -0.011 0.011 0.011 -0.011 -0.004 0.004 -0.000 0.004 0.004 -0.011 -0.011 0.011 0.011 -0.004 -0.004 1 2 3 4 5 6 7 8 9 0.041 0.004 -0.004 -0.011 0.011 0.011 -0.011 -0.004 0.004 -0.000 -0.004 -0.004 0.011 0.011 -0.011 -0.011 0.004 0.004 Upper and lower wings at the 2 th omega 0.0000 16.7000 [eV] 1 2 3 4 5 6 7 8 9 0.484 0.007 -0.007 -0.022 0.022 0.022 -0.022 -0.008 0.008 -0.000 0.007 0.007 -0.022 -0.022 0.022 0.022 -0.008 -0.008 1 2 3 4 5 6 7 8 9 0.484 0.007 -0.007 -0.022 0.022 0.022 -0.022 -0.008 0.008 -0.000 -0.007 -0.007 0.022 0.022 -0.022 -0.022 0.008 0.008 -------------------------------------------------------------------------------- q-point number 2 q = ( 0.500000, 0.500000, 0.000000) [r.l.u.] -------------------------------------------------------------------------------- chi0(G,G') at the 1 th omega 0.0000 0.0000 [eV] 1 2 3 4 5 6 7 8 9 1 -18.011 -3.066 -1.214 -1.217 -3.062 -1.217 -3.062 -3.063 -1.210 -0.000 -3.066 1.214 -1.217 3.062 -1.217 3.062 -3.063 1.210 2 -3.066 -12.569 0.000 -0.019 -0.000 -0.029 0.000 -0.067 -0.000 3.066 -0.000 -2.507 -0.000 0.598 0.000 0.610 -0.000 -0.226 chi0(G,G') at the 2 th omega 0.0000 16.7000 [eV] 1 2 3 4 5 6 7 8 9 1 -4.167 -0.794 -0.285 -0.287 -0.792 -0.287 -0.792 -0.792 -0.282 -0.000 -0.794 0.285 -0.287 0.792 -0.287 0.792 -0.792 0.282 2 -0.794 -4.302 0.000 -0.042 -0.000 -0.048 0.000 0.055 -0.000 0.794 0.000 -0.577 -0.000 -0.002 0.000 0.007 -0.000 -0.004 Average fulfillment of the sum rule on Im[epsilon] for q-point 2 : 59.28 [%] -------------------------------------------------------------------------------- q-point number 3 q = ( 0.500000, 0.250000, 0.250000) [r.l.u.] -------------------------------------------------------------------------------- chi0(G,G') at the 1 th omega 0.0000 0.0000 [eV] 1 2 3 4 5 6 7 8 9 1 -14.681 -2.448 0.459 -2.350 -2.351 0.459 -2.448 -2.354 -2.352 -0.000 -2.448 -0.459 -2.350 2.351 0.459 2.448 -2.354 2.352 2 -2.448 -12.885 0.000 -0.422 0.000 -0.028 -0.000 -0.424 -0.000 2.448 0.000 -3.364 0.000 0.186 -0.000 0.460 -0.000 0.171 chi0(G,G') at the 2 th omega 0.0000 16.7000 [eV] 1 2 3 4 5 6 7 8 9 1 -2.545 -0.657 0.234 -0.426 -0.427 0.234 -0.657 -0.428 -0.427 0.000 -0.657 -0.234 -0.426 0.427 0.234 0.657 -0.428 0.427 2 -0.657 -4.545 0.000 0.022 0.000 0.011 -0.000 0.019 -0.000 0.657 0.000 -0.825 0.000 0.005 -0.000 0.047 -0.000 -0.003 Average fulfillment of the sum rule on Im[epsilon] for q-point 3 : 63.90 [%] -------------------------------------------------------------------------------- q-point number 4 q = ( 0.000000, 0.500000, 0.000000) [r.l.u.] -------------------------------------------------------------------------------- chi0(G,G') at the 1 th omega 0.0000 0.0000 [eV] 1 2 3 4 5 6 7 8 9 1 -14.476 -2.047 -2.646 -2.046 -2.645 -2.046 -2.645 -2.108 3.151 -0.000 -2.047 2.646 -2.046 2.645 -2.046 2.645 -2.108 -3.151 2 -2.047 -17.857 0.000 0.354 -0.000 0.338 0.000 0.264 -0.000 2.047 -0.000 -2.604 -0.000 0.245 0.000 0.247 0.000 -0.918 chi0(G,G') at the 2 th omega 0.0000 16.7000 [eV] 1 2 3 4 5 6 7 8 9 1 -3.449 -0.388 -0.687 -0.387 -0.687 -0.387 -0.687 -0.780 0.524 -0.000 -0.388 0.687 -0.387 0.687 -0.387 0.687 -0.780 -0.524 2 -0.388 -5.661 0.000 -0.021 -0.000 -0.032 0.000 0.038 -0.000 0.388 -0.000 -0.611 -0.000 -0.027 0.000 -0.025 0.000 -0.112 Average fulfillment of the sum rule on Im[epsilon] for q-point 4 : 62.03 [%] -------------------------------------------------------------------------------- q-point number 5 q = ( 0.500000,-0.250000, 0.250000) [r.l.u.] -------------------------------------------------------------------------------- chi0(G,G') at the 1 th omega 0.0000 0.0000 [eV] 1 2 3 4 5 6 7 8 9 1 -18.994 -2.734 -1.405 -2.887 -2.489 -2.483 -2.895 -1.406 -2.735 -0.000 -2.734 1.405 -2.887 2.489 -2.483 2.895 -1.406 2.735 2 -2.734 -10.375 0.000 0.266 -0.000 0.486 0.000 -0.353 -0.000 2.734 0.000 -2.434 -0.000 0.305 0.000 0.409 -0.000 0.167 chi0(G,G') at the 2 th omega 0.0000 16.7000 [eV] 1 2 3 4 5 6 7 8 9 1 -4.682 -0.838 -0.093 -0.808 -0.665 -0.662 -0.813 -0.093 -0.839 -0.000 -0.838 0.093 -0.808 0.665 -0.662 0.813 -0.093 0.839 2 -0.838 -3.586 0.000 0.094 -0.000 0.040 0.000 -0.023 -0.000 0.838 -0.000 -0.492 -0.000 0.034 0.000 0.033 0.000 0.017 Average fulfillment of the sum rule on Im[epsilon] for q-point 5 : 59.51 [%] -------------------------------------------------------------------------------- q-point number 6 q = ( 0.000000,-0.250000,-0.250000) [r.l.u.] -------------------------------------------------------------------------------- chi0(G,G') at the 1 th omega 0.0000 0.0000 [eV] 1 2 3 4 5 6 7 8 9 1 -10.807 -2.286 -0.057 -0.056 -2.285 -2.285 -0.057 -0.057 -2.286 0.000 -2.286 0.057 -0.056 2.285 -2.285 0.057 -0.057 2.286 2 -2.286 -14.835 -0.000 -0.406 0.000 -0.454 0.000 -0.435 -0.000 2.286 0.000 -3.485 0.000 0.281 -0.000 -0.365 0.000 0.277 chi0(G,G') at the 2 th omega 0.0000 16.7000 [eV] 1 2 3 4 5 6 7 8 9 1 -1.412 -0.428 0.076 0.076 -0.428 -0.428 0.076 0.075 -0.429 -0.000 -0.428 -0.076 0.076 0.428 -0.428 -0.076 0.075 0.429 2 -0.428 -5.190 -0.000 0.031 0.000 0.036 0.000 0.011 -0.000 0.428 0.000 -0.910 0.000 0.016 -0.000 -0.090 0.000 0.013 Average fulfillment of the sum rule on Im[epsilon] for q-point 6 : 71.01 [%] ================================================================================ == DATASET 4 ================================================================== - mpi_nproc: 1, omp_nthreads: -1 (-1 if OMP is not activated) --- !DatasetInfo iteration_state: {dtset: 4, } dimensions: {natom: 2, nkpt: 6, mband: 30, nsppol: 1, nspinor: 1, nspden: 1, mpw: 266, } cutoff_energies: {ecut: 7.6, pawecutdg: -1.0, } electrons: {nelect: 8.00000000E+00, charge: 0.00000000E+00, occopt: 1.00000000E+00, tsmear: 1.00000000E-02, } meta: {optdriver: 4, gwcalctyp: 0, } ... mkfilename : getwfk/=0, take file _WFK from output of DATASET 2. mkfilename : getscr/=0, take file _SCR from output of DATASET 3. Exchange-correlation functional for the present dataset will be: LDA: new Teter (4/93) with spin-polarized option - ixc=1 Citation for XC functional: S. Goedecker, M. Teter, J. Huetter, PRB 54, 1703 (1996) SIGMA: Calculation of the GW corrections Based on a program developped by R.W. Godby, V. Olevano, G. Onida, and L. Reining. Incorporated in ABINIT by V. Olevano, G.-M. Rignanese, and M. Torrent. .Using double precision arithmetic ; gwpc = 8 Real(R)+Recip(G) space primitive vectors, cartesian coordinates (Bohr,Bohr^-1): R(1)= 0.0000000 5.1085000 5.1085000 G(1)= -0.0978761 0.0978761 0.0978761 R(2)= 5.1085000 0.0000000 5.1085000 G(2)= 0.0978761 -0.0978761 0.0978761 R(3)= 5.1085000 5.1085000 0.0000000 G(3)= 0.0978761 0.0978761 -0.0978761 Unit cell volume ucvol= 2.6663072E+02 bohr^3 Angles (23,13,12)= 6.00000000E+01 6.00000000E+01 6.00000000E+01 degrees -------------------------------------------------------------------------------- ==== K-mesh for the wavefunctions ==== Number of points in the irreducible wedge : 6 Reduced coordinates and weights : 1) -2.50000000E-01 -2.50000000E-01 0.00000000E+00 0.18750 2) -2.50000000E-01 2.50000000E-01 0.00000000E+00 0.37500 3) 5.00000000E-01 5.00000000E-01 0.00000000E+00 0.09375 4) -2.50000000E-01 5.00000000E-01 2.50000000E-01 0.18750 5) 5.00000000E-01 0.00000000E+00 0.00000000E+00 0.12500 6) 0.00000000E+00 0.00000000E+00 0.00000000E+00 0.03125 Together with 48 symmetry operations and time-reversal symmetry yields 32 points in the full Brillouin Zone. ==== Q-mesh for screening function ==== Number of points in the irreducible wedge : 6 Reduced coordinates and weights : 1) 0.00000000E+00 0.00000000E+00 0.00000000E+00 0.03125 2) 5.00000000E-01 5.00000000E-01 0.00000000E+00 0.09375 3) 5.00000000E-01 2.50000000E-01 2.50000000E-01 0.37500 4) 0.00000000E+00 5.00000000E-01 0.00000000E+00 0.12500 5) 5.00000000E-01 -2.50000000E-01 2.50000000E-01 0.18750 6) 0.00000000E+00 -2.50000000E-01 -2.50000000E-01 0.18750 Together with 48 symmetry operations and time-reversal symmetry yields 32 points in the full Brillouin Zone. setmesh: FFT mesh size selected = 20x 20x 20 total number of points = 8000 Number of electrons calculated from density = 7.9999; Expected = 8.0000 average of density, n = 0.030004 r_s = 1.9964 omega_plasma = 16.7087 [eV] === KS Band Gaps === >>>> For spin 1 Minimum direct gap = 2.5053 [eV], located at k-point : 0.0000 0.0000 0.0000 Fundamental gap = 0.5581 [eV], Top of valence bands at : 0.0000 0.0000 0.0000 Bottom of conduction at : 0.5000 0.5000 0.0000 SIGMA fundamental parameters: PLASMON POLE MODEL 1 number of plane-waves for SigmaX 283 number of plane-waves for SigmaC and W 89 number of plane-waves for wavefunctions 283 number of bands 30 number of independent spin polarizations 1 number of spinorial components 1 number of k-points in IBZ 6 number of q-points in IBZ 6 number of symmetry operations 48 number of k-points in BZ 32 number of q-points in BZ 32 number of frequencies for dSigma/dE 9 frequency step for dSigma/dE [eV] 0.25 number of omega for Sigma on real axis 0 max omega for Sigma on real axis [eV] 0.00 zcut for avoiding poles [eV] 0.10 EPSILON^-1 parameters (SCR file): dimension of the eps^-1 matrix on file 89 dimension of the eps^-1 matrix used 89 number of plane-waves for wavefunctions 283 number of bands 17 number of q-points in IBZ 6 number of frequencies 2 number of real frequencies 1 number of imag frequencies 1 matrix elements of self-energy operator (all in [eV]) Perturbative Calculation --- !SelfEnergy_ee iteration_state: {dtset: 4, } kpoint : [ 0.000, 0.000, 0.000, ] spin : 1 KS_gap : 2.505 QP_gap : 3.119 Delta_QP_KS: 0.614 data: !SigmaeeData | Band E0 <VxcDFT> SigX SigC(E0) Z dSigC/dE Sig(E) E-E0 E 2 5.967 -11.268 -13.253 1.814 0.770 -0.299 -11.399 -0.131 5.836 3 5.967 -11.268 -13.253 1.814 0.770 -0.299 -11.399 -0.131 5.836 4 5.967 -11.268 -13.253 1.814 0.770 -0.299 -11.399 -0.131 5.836 5 8.472 -10.056 -5.573 -3.856 0.771 -0.298 -9.573 0.483 8.955 6 8.472 -10.056 -5.573 -3.856 0.771 -0.298 -9.573 0.483 8.955 7 8.472 -10.056 -5.573 -3.856 0.771 -0.298 -9.573 0.483 8.955 ... == END DATASET(S) ============================================================== ================================================================================ -outvars: echo values of variables after computation -------- acell 1.0217000000E+01 1.0217000000E+01 1.0217000000E+01 Bohr amu 2.80855000E+01 bdgw4 4 5 diemac 1.20000000E+01 ecut1 8.00000000E+00 Hartree ecut2 8.00000000E+00 Hartree ecut3 7.56385066E+00 Hartree ecut4 7.56385066E+00 Hartree ecuteps1 0.00000000E+00 Hartree ecuteps2 0.00000000E+00 Hartree ecuteps3 3.59282906E+00 Hartree ecuteps4 0.00000000E+00 Hartree ecutsigx1 0.00000000E+00 Hartree ecutsigx2 0.00000000E+00 Hartree ecutsigx3 0.00000000E+00 Hartree ecutsigx4 7.56385066E+00 Hartree ecutwfn1 0.00000000E+00 Hartree ecutwfn2 0.00000000E+00 Hartree ecutwfn3 7.56385066E+00 Hartree ecutwfn4 7.56385066E+00 Hartree etotal1 -8.8642542742E+00 etotal3 0.0000000000E+00 etotal4 0.0000000000E+00 fcart1 -0.0000000000E+00 -0.0000000000E+00 -0.0000000000E+00 -0.0000000000E+00 -0.0000000000E+00 -0.0000000000E+00 fcart3 0.0000000000E+00 0.0000000000E+00 0.0000000000E+00 0.0000000000E+00 0.0000000000E+00 0.0000000000E+00 fcart4 0.0000000000E+00 0.0000000000E+00 0.0000000000E+00 0.0000000000E+00 0.0000000000E+00 0.0000000000E+00 - fftalg 312 getden1 0 getden2 -1 getden3 0 getden4 0 getscr1 0 getscr2 0 getscr3 0 getscr4 -1 getwfk1 0 getwfk2 0 getwfk3 -1 getwfk4 -2 iscf1 7 iscf2 -2 iscf3 7 iscf4 7 istwfk 0 0 1 0 1 1 jdtset 1 2 3 4 kpt -2.50000000E-01 -2.50000000E-01 0.00000000E+00 -2.50000000E-01 2.50000000E-01 0.00000000E+00 5.00000000E-01 5.00000000E-01 0.00000000E+00 -2.50000000E-01 5.00000000E-01 2.50000000E-01 5.00000000E-01 0.00000000E+00 0.00000000E+00 0.00000000E+00 0.00000000E+00 0.00000000E+00 kptgw4 0.00000000E+00 0.00000000E+00 0.00000000E+00 kptrlatt 2 -2 2 -2 2 2 -2 -2 2 kptrlen 2.04340000E+01 P mkmem 6 natom 2 nband1 6 nband2 40 nband3 17 nband4 30 nbdbuf1 0 nbdbuf2 5 nbdbuf3 0 nbdbuf4 0 ndtset 4 ngfft1 20 20 20 ngfft2 20 20 20 ngfft3 18 18 18 ngfft4 18 18 18 nkpt 6 nkptgw1 0 nkptgw2 0 nkptgw3 0 nkptgw4 1 npweps1 0 npweps2 0 npweps3 89 npweps4 0 npwsigx1 0 npwsigx2 0 npwsigx3 0 npwsigx4 283 npwwfn1 0 npwwfn2 0 npwwfn3 283 npwwfn4 283 nstep 20 nsym 48 ntypat 1 occ1 2.000000 2.000000 2.000000 2.000000 0.000000 0.000000 occ3 2.000000 2.000000 2.000000 2.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 occ4 2.000000 2.000000 2.000000 2.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 optdriver1 0 optdriver2 0 optdriver3 3 optdriver4 4 ppmfrq1 0.00000000E+00 Hartree ppmfrq2 0.00000000E+00 Hartree ppmfrq3 6.13713734E-01 Hartree ppmfrq4 0.00000000E+00 Hartree rprim 0.0000000000E+00 5.0000000000E-01 5.0000000000E-01 5.0000000000E-01 0.0000000000E+00 5.0000000000E-01 5.0000000000E-01 5.0000000000E-01 0.0000000000E+00 spgroup 227 strten1 9.8731950055E-07 9.8731950055E-07 9.8731950055E-07 0.0000000000E+00 0.0000000000E+00 0.0000000000E+00 strten3 0.0000000000E+00 0.0000000000E+00 0.0000000000E+00 0.0000000000E+00 0.0000000000E+00 0.0000000000E+00 strten4 0.0000000000E+00 0.0000000000E+00 0.0000000000E+00 0.0000000000E+00 0.0000000000E+00 0.0000000000E+00 symrel 1 0 0 0 1 0 0 0 1 -1 0 0 0 -1 0 0 0 -1 0 -1 1 0 -1 0 1 -1 0 0 1 -1 0 1 0 -1 1 0 -1 0 0 -1 0 1 -1 1 0 1 0 0 1 0 -1 1 -1 0 0 1 -1 1 0 -1 0 0 -1 0 -1 1 -1 0 1 0 0 1 -1 0 0 -1 1 0 -1 0 1 1 0 0 1 -1 0 1 0 -1 0 -1 1 1 -1 0 0 -1 0 0 1 -1 -1 1 0 0 1 0 1 0 0 0 0 1 0 1 0 -1 0 0 0 0 -1 0 -1 0 0 1 -1 0 0 -1 1 0 -1 0 -1 1 0 0 1 -1 0 1 -1 0 1 -1 1 0 -1 0 0 1 0 -1 1 -1 0 1 0 0 0 -1 0 1 -1 0 0 -1 1 0 1 0 -1 1 0 0 1 -1 1 0 -1 0 0 -1 0 1 -1 -1 0 1 0 0 1 0 -1 1 0 1 0 0 0 1 1 0 0 0 -1 0 0 0 -1 -1 0 0 1 0 -1 0 1 -1 0 0 -1 -1 0 1 0 -1 1 0 0 1 0 -1 0 0 -1 1 1 -1 0 0 1 0 0 1 -1 -1 1 0 -1 0 1 -1 0 0 -1 1 0 1 0 -1 1 0 0 1 -1 0 0 1 0 1 0 0 0 0 1 0 -1 0 -1 0 0 0 0 -1 0 0 -1 0 1 -1 1 0 -1 0 0 1 0 -1 1 -1 0 1 1 -1 0 0 -1 1 0 -1 0 -1 1 0 0 1 -1 0 1 0 0 0 1 1 0 0 0 1 0 0 0 -1 -1 0 0 0 -1 0 -1 1 0 -1 0 0 -1 0 1 1 -1 0 1 0 0 1 0 -1 0 0 1 0 1 0 1 0 0 0 0 -1 0 -1 0 -1 0 0 1 -1 0 0 -1 0 0 -1 1 -1 1 0 0 1 0 0 1 -1 0 0 -1 1 0 -1 0 1 -1 0 0 1 -1 0 1 0 -1 1 -1 1 0 -1 0 1 -1 0 0 1 -1 0 1 0 -1 1 0 0 tnons 0.0000000 0.0000000 0.0000000 0.2500000 0.2500000 0.2500000 0.0000000 0.0000000 0.0000000 0.2500000 0.2500000 0.2500000 0.0000000 0.0000000 0.0000000 0.2500000 0.2500000 0.2500000 0.0000000 0.0000000 0.0000000 0.2500000 0.2500000 0.2500000 0.0000000 0.0000000 0.0000000 0.2500000 0.2500000 0.2500000 0.0000000 0.0000000 0.0000000 0.2500000 0.2500000 0.2500000 0.0000000 0.0000000 0.0000000 0.2500000 0.2500000 0.2500000 0.0000000 0.0000000 0.0000000 0.2500000 0.2500000 0.2500000 0.0000000 0.0000000 0.0000000 0.2500000 0.2500000 0.2500000 0.0000000 0.0000000 0.0000000 0.2500000 0.2500000 0.2500000 0.0000000 0.0000000 0.0000000 0.2500000 0.2500000 0.2500000 0.0000000 0.0000000 0.0000000 0.2500000 0.2500000 0.2500000 0.0000000 0.0000000 0.0000000 0.2500000 0.2500000 0.2500000 0.0000000 0.0000000 0.0000000 0.2500000 0.2500000 0.2500000 0.0000000 0.0000000 0.0000000 0.2500000 0.2500000 0.2500000 0.0000000 0.0000000 0.0000000 0.2500000 0.2500000 0.2500000 0.0000000 0.0000000 0.0000000 0.2500000 0.2500000 0.2500000 0.0000000 0.0000000 0.0000000 0.2500000 0.2500000 0.2500000 0.0000000 0.0000000 0.0000000 0.2500000 0.2500000 0.2500000 0.0000000 0.0000000 0.0000000 0.2500000 0.2500000 0.2500000 0.0000000 0.0000000 0.0000000 0.2500000 0.2500000 0.2500000 0.0000000 0.0000000 0.0000000 0.2500000 0.2500000 0.2500000 0.0000000 0.0000000 0.0000000 0.2500000 0.2500000 0.2500000 0.0000000 0.0000000 0.0000000 0.2500000 0.2500000 0.2500000 tolvrs1 1.00000000E-10 tolvrs2 0.00000000E+00 tolvrs3 0.00000000E+00 tolvrs4 0.00000000E+00 tolwfr1 0.00000000E+00 tolwfr2 1.00000000E-18 tolwfr3 0.00000000E+00 tolwfr4 0.00000000E+00 typat 1 1 wtk 0.18750 0.37500 0.09375 0.18750 0.12500 0.03125 xangst 0.0000000000E+00 0.0000000000E+00 0.0000000000E+00 1.3516508850E+00 1.3516508850E+00 1.3516508850E+00 xcart 0.0000000000E+00 0.0000000000E+00 0.0000000000E+00 2.5542500000E+00 2.5542500000E+00 2.5542500000E+00 xred 0.0000000000E+00 0.0000000000E+00 0.0000000000E+00 2.5000000000E-01 2.5000000000E-01 2.5000000000E-01 znucl 14.00000 ================================================================================ - Timing analysis has been suppressed with timopt=0 ================================================================================ Suggested references for the acknowledgment of ABINIT usage. The users of ABINIT have little formal obligations with respect to the ABINIT group (those specified in the GNU General Public License, http://www.gnu.org/copyleft/gpl.txt). However, it is common practice in the scientific literature, to acknowledge the efforts of people that have made the research possible. In this spirit, please find below suggested citations of work written by ABINIT developers, corresponding to implementations inside of ABINIT that you have used in the present run. Note also that it will be of great value to readers of publications presenting these results, to read papers enabling them to understand the theoretical formalism and details of the ABINIT implementation. For information on why they are suggested, see also https://docs.abinit.org/theory/acknowledgments. - - [1] The Abinit project: Impact, environment and recent developments. - Computer Phys. Comm. 248, 107042 (2020). - X.Gonze, B. Amadon, G. Antonius, F.Arnardi, L.Baguet, J.-M.Beuken, - J.Bieder, F.Bottin, J.Bouchet, E.Bousquet, N.Brouwer, F.Bruneval, - G.Brunin, T.Cavignac, J.-B. Charraud, Wei Chen, M.Cote, S.Cottenier, - J.Denier, G.Geneste, Ph.Ghosez, M.Giantomassi, Y.Gillet, O.Gingras, - D.R.Hamann, G.Hautier, Xu He, N.Helbig, N.Holzwarth, Y.Jia, F.Jollet, - W.Lafargue-Dit-Hauret, K.Lejaeghere, M.A.L.Marques, A.Martin, C.Martins, - H.P.C. Miranda, F.Naccarato, K. Persson, G.Petretto, V.Planes, Y.Pouillon, - S.Prokhorenko, F.Ricci, G.-M.Rignanese, A.H.Romero, M.M.Schmitt, M.Torrent, - M.J.van Setten, B.Van Troeye, M.J.Verstraete, G.Zerah and J.W.Zwanzig - Comment: the fifth generic paper describing the ABINIT project. - Note that a version of this paper, that is not formatted for Computer Phys. Comm. - is available at https://www.abinit.org/sites/default/files/ABINIT20.pdf . - The licence allows the authors to put it on the Web. - DOI and bibtex: see https://docs.abinit.org/theory/bibliography/#gonze2020 - - [2] Recent developments in the ABINIT software package. - Computer Phys. Comm. 205, 106 (2016). - X.Gonze, F.Jollet, F.Abreu Araujo, D.Adams, B.Amadon, T.Applencourt, - C.Audouze, J.-M.Beuken, J.Bieder, A.Bokhanchuk, E.Bousquet, F.Bruneval - D.Caliste, M.Cote, F.Dahm, F.Da Pieve, M.Delaveau, M.Di Gennaro, - B.Dorado, C.Espejo, G.Geneste, L.Genovese, A.Gerossier, M.Giantomassi, - Y.Gillet, D.R.Hamann, L.He, G.Jomard, J.Laflamme Janssen, S.Le Roux, - A.Levitt, A.Lherbier, F.Liu, I.Lukacevic, A.Martin, C.Martins, - M.J.T.Oliveira, S.Ponce, Y.Pouillon, T.Rangel, G.-M.Rignanese, - A.H.Romero, B.Rousseau, O.Rubel, A.A.Shukri, M.Stankovski, M.Torrent, - M.J.Van Setten, B.Van Troeye, M.J.Verstraete, D.Waroquier, J.Wiktor, - B.Xu, A.Zhou, J.W.Zwanziger. - Comment: the fourth generic paper describing the ABINIT project. - Note that a version of this paper, that is not formatted for Computer Phys. Comm. - is available at https://www.abinit.org/sites/default/files/ABINIT16.pdf . - The licence allows the authors to put it on the Web. - DOI and bibtex: see https://docs.abinit.org/theory/bibliography/#gonze2016 - - [3] ABINIT: First-principles approach of materials and nanosystem properties. - Computer Phys. Comm. 180, 2582-2615 (2009). - X. Gonze, B. Amadon, P.-M. Anglade, J.-M. Beuken, F. Bottin, P. Boulanger, F. Bruneval, - D. Caliste, R. Caracas, M. Cote, T. Deutsch, L. Genovese, Ph. Ghosez, M. Giantomassi - S. Goedecker, D.R. Hamann, P. Hermet, F. Jollet, G. Jomard, S. Leroux, M. Mancini, S. Mazevet, - M.J.T. Oliveira, G. Onida, Y. Pouillon, T. Rangel, G.-M. Rignanese, D. Sangalli, R. Shaltaf, - M. Torrent, M.J. Verstraete, G. Zerah, J.W. Zwanziger - Comment: the third generic paper describing the ABINIT project. - Note that a version of this paper, that is not formatted for Computer Phys. Comm. - is available at https://www.abinit.org/sites/default/files/ABINIT_CPC_v10.pdf . - The licence allows the authors to put it on the Web. - DOI and bibtex: see https://docs.abinit.org/theory/bibliography/#gonze2009 - - And optionally: - - [4] A brief introduction to the ABINIT software package. - Z. Kristallogr. 220, 558-562 (2005). - X. Gonze, G.-M. Rignanese, M. Verstraete, J.-M. Beuken, Y. Pouillon, R. Caracas, F. Jollet, - M. Torrent, G. Zerah, M. Mikami, Ph. Ghosez, M. Veithen, J.-Y. Raty, V. Olevano, F. Bruneval, - L. Reining, R. Godby, G. Onida, D.R. Hamann, and D.C. Allan. - Comment: the second generic paper describing the ABINIT project. Note that this paper - should be cited especially if you are using the GW part of ABINIT, as several authors - of this part are not in the list of authors of the first or third paper. - The .pdf of the latter paper is available at https://www.abinit.org/sites/default/files/zfk_0505-06_558-562.pdf. - Note that it should not redistributed (Copyright by Oldenburg Wissenschaftverlag, - the licence allows the authors to put it on the Web). - DOI and bibtex: see https://docs.abinit.org/theory/bibliography/#gonze2005 - - Proc. 0 individual time (sec): cpu= 6.8 wall= 17.5 ================================================================================ Calculation completed. .Delivered 5 WARNINGs and 7 COMMENTs to log file. +Overall time at end (sec) : cpu= 6.8 wall= 17.5
After the description of the unit cell and of the pseudopotentials, you will find the list of k-points used for the electrons and the grid of q-points (in the Irreducible part of the Brillouin Zone) on which the susceptibility and dielectric matrices will be computed.
==== K-mesh for the wavefunctions ==== Number of points in the irreducible wedge : 6 Reduced coordinates and weights : 1) -2.50000000E-01 -2.50000000E-01 0.00000000E+00 0.18750 2) -2.50000000E-01 2.50000000E-01 0.00000000E+00 0.37500 3) 5.00000000E-01 5.00000000E-01 0.00000000E+00 0.09375 4) -2.50000000E-01 5.00000000E-01 2.50000000E-01 0.18750 5) 5.00000000E-01 0.00000000E+00 0.00000000E+00 0.12500 6) 0.00000000E+00 0.00000000E+00 0.00000000E+00 0.03125 Together with 48 symmetry operations and time-reversal symmetry yields 32 points in the full Brillouin Zone. ==== Q-mesh for the screening function ==== Number of points in the irreducible wedge : 6 Reduced coordinates and weights : 1) 0.00000000E+00 0.00000000E+00 0.00000000E+00 0.03125 2) 5.00000000E-01 5.00000000E-01 0.00000000E+00 0.09375 3) 5.00000000E-01 2.50000000E-01 2.50000000E-01 0.37500 4) 0.00000000E+00 5.00000000E-01 0.00000000E+00 0.12500 5) 5.00000000E-01 -2.50000000E-01 2.50000000E-01 0.18750 6) 0.00000000E+00 -2.50000000E-01 -2.50000000E-01 0.18750 Together with 48 symmetry operations and time-reversal symmetry yields 32 points in the full Brillouin Zone.
The q-mesh is the set of points defined as all the possible differences among the k-points ( \mathbf{q} =\mathbf{k}-\mathbf{k}' ) of the grid chosen to generate the WFK file. From the last statement it is clear the importance of choosing homogeneous k-point grids in order to minimize the number of q-points is clear.
After this section, the code prints the parameters of the FFT grid needed to represent the wavefunctions and to compute their convolution (required for the screening matrices). Then we have some information about the MPI distribution of the bands and the total number of valence electrons computed by integrating the density in the unit cell.
setmesh: FFT mesh size selected = 20x 20x 20 total number of points = 8000 - screening: taking advantage of time-reversal symmetry - Maximum band index for partially occupied states nbvw = 4 - Remaining bands to be divided among processors nbcw = 13 - Number of bands treated by each node ~13 Number of electrons calculated from density = 7.9999; Expected = 8.0000 average of density, n = 0.030004 r_s = 1.9964 omega_plasma = 16.7087 [eV]
On the basis of the density, one can obtain the classical Drude plasmon frequency. The next lines calculate the average density of the system, and evaluate the Wigner radius r_s, then compute the Drude plasmon frequency.
Number of electrons calculated from density = 7.9999; Expected = 8.0000 average of density, n = 0.030004 r_s = 1.9964 omega_plasma = 16.7087 [eV]
This is the value used by default for ppmfrq. It is in fact the second frequency where the code calculates the dielectric matrix to adjust the plasmon-pole model parameters.
It has been found that Drude plasma frequency is a reasonable value where to adjust the model. The control over this parameter is however left to the user in order to check that the result does not change when changing ppmfrq. If it is the case, then the plasmon-pole model is not appropriated and one should go beyond by taking into account a full dynamical dependence in the screening (see later, the contour-deformation method). However, the plasmon-pole model has been found to work well for a very large range of solid-state systems when focusing only on the real part of the GW corrections.
At the end of the screening calculation, the macroscopic dielectric constant is printed:
dielectric constant = 24.4980 dielectric constant without local fields = 27.1361
Note
Note that the convergence in the dielectric constant does not guarantee the convergence in the GW corrections. In fact, the dielectric constant is representative of only one element i.e. the head of the full dielectric matrix. Even if the convergence on the dielectric constant with local fields takes somehow into account also other non-diagonal elements. In a GW calculation the whole \epsilon^{-1} matrix is used to build the Self-Energy operator.
The dielectric constant reported here is the so-called RPA dielectric constant due to the electrons. Although evaluated at zero frequency, it is understood that the ionic response is not included (this term can be computed with DFPT and ANADDB). The RPA dielectric constant restricted to electronic effects is also not the same as the one computed in the DFPT part of ABINIT, that includes exchange-correlation effects.
We now enter the fourth dataset. As for dataset 3, after some general information (origin of WFK file, header, description of unit cell, k-points, q-points), the description of the FFT grid and jellium parameters, there is the echo of parameters for the plasmon-pole model, and the inverse dielectric function (the screening). The self-energy operator has been constructed, and one can evaluate the GW energies for each state.
The final results are:
--- !SelfEnergy_ee kpoint : [ 0.000, 0.000, 0.000, ] spin : 1 KS_gap : 2.505 QP_gap : 3.120 Delta_QP_KS: 0.614 data: !SigmaeeData | Band E0 <VxcDFT> SigX SigC(E0) Z dSigC/dE Sig(E) E-E0 E 4 5.967 -11.268 -13.253 1.814 0.770 -0.299 -11.400 -0.132 5.835 5 8.472 -10.056 -5.573 -3.856 0.770 -0.298 -9.573 0.483 8.955 ...
For the desired k-point (\Gamma point), for state 4, then state 5, one finds different information:
- E0 is the KS eigenenergy
- VxcDFT gives the average KS exchange-correlation potential
- SigX gives the exchange contribution to the self-energy
- SigC(E0) gives the correlation contribution to the self-energy, evaluated at the KS eigenenergy
- Z is the renormalisation factor
- dSigC/dE is the energy derivative of SigC with respect to the energy
- SigC(E) gives the correlation contribution to the self-energy, evaluated at the GW energy
- E-E0 is the difference between GW energy and KS eigenenergy
- E is the GW quasiparticle energy
In this case, the gap is also analyzed: E^0_gap is the direct KS gap at that particular k-point (and spin, in the case of spin-polarized calculations), E^GW_gap is the GW one, and DeltaE^GW_gap is the difference. This direct gap is always computed between the band whose number is equal to the number of electrons in the cell divided by two (integer part, in case of spin-polarized calculation), and the next one. This means that the value reported by the code may be wrong if the final QP energies obtained in the perturbative approach are not ordered by increasing energy anymore. So it’s always a good idea to check that the “gap” reported by the code corresponds to the real QP direct gap.
Warning
For a metal, these two bands do not systematically lie below and above the KS Fermi energy - but the concept of a direct gap is not relevant in that case. Moreover one should compute the Fermy energy of the QP system.
It is seen that the average KS exchange-correlation potential for the state 4 (a valence state) is rather close to the exchange self-energy correction. For that state, the correlation correction is small, and the difference between KS and GW energies is also small (0.128 eV). By contrast, the exchange self-energy is much smaller than the average Kohn-Sham potential for the state 5 (a conduction state), but the correlation correction is much larger than for state 4. On the whole, the difference between Kohn- Sham and GW energies is not very large, but nevertheless, it is quite important when compared with the size of the gap.
If AbiPy is installed on your machine, you can use the abiopen script
with the --print
option to extract the results from the SIGRES.nc file
and print them to terminal:
abiopen.py tgw1_1o_DS4_SIGRES.nc -p ================================= Structure ================================= Full Formula (Si2) Reduced Formula: Si abc : 3.823046 3.823046 3.823046 angles: 60.000000 60.000000 60.000000 Sites (2) # SP a b c --- ---- ---- ---- ---- 0 Si 0 0 0 1 Si 0.25 0.25 0.25 Abinit Spacegroup: spgid: 0, num_spatial_symmetries: 48, has_timerev: True, symmorphic: True ============================== Kohn-Sham bands ============================== Number of electrons: 8.0, Fermi level: 6.246 (eV) nsppol: 1, nkpt: 6, mband: 30, nspinor: 1, nspden: 1 smearing scheme: none, tsmear_eV: 0.272, occopt: 1 Direct gap: Energy: 2.505 (eV) Initial state: spin=0, kpt=[+0.000, +0.000, +0.000], weight: 0.031, band=3, eig=5.967, occ=2.000 Final state: spin=0, kpt=[+0.000, +0.000, +0.000], weight: 0.031, band=4, eig=8.472, occ=0.000 Fundamental gap: Energy: 0.558 (eV) Initial state: spin=0, kpt=[+0.000, +0.000, +0.000], weight: 0.031, band=3, eig=5.967, occ=2.000 Final state: spin=0, kpt=[+0.500, +0.500, +0.000], weight: 0.094, band=4, eig=6.525, occ=0.000 Bandwidth: 12.101 (eV) Valence maximum located at: spin=0, kpt=[+0.000, +0.000, +0.000], weight: 0.031, band=3, eig=5.967, occ=2.000 Conduction minimum located at: spin=0, kpt=[+0.500, +0.500, +0.000], weight: 0.094, band=4, eig=6.525, occ=0.000 =============================== QP direct gaps =============================== QP_dirgap: 3.120 (eV) for K-point: [+0.000, +0.000, +0.000] $\Gamma$, spin: 0 ============== QP results for each k-point and spin (All in eV) ============== K-point: [+0.000, +0.000, +0.000] $\Gamma$, spin: 0 band e0 qpe qpe_diago vxcme sigxme sigcmee0 vUme ze0 3 3 5.967 5.835 5.796 -11.268 -13.253 1.814 0.0 0.77 4 4 8.472 8.955 9.099 -10.056 -5.573 -3.856 0.0 0.77
For further details about the SIGRES.nc file and the AbiPy API see the Sigres notebook .
2 Preparing convergence studies: Kohn-Sham structure (WFK file) and screening (SCR file)¶
In the following sections, we will perform different convergence studies. In order to keep the CPU time at a reasonable level, we will use fake WFK and SCR data. Moreover we will only consider the correction at the \Gamma point only. In this way, we will be able to verify convergence aspects that could be very cumbersome (at least in the framework of a tutorial) if more k-points were used. Testing the convergence with a \Gamma point only grid of k-point represents a convenient approach although some caution should always be used.
In directory Work_gw1, copy the file ../tgw1_2.in, and modify the tgw1_x.files file as usual. Edit the tgw1_2.in file, and take the time to examine it.
Then, issue:
abinit < tgw1_x.files > tgw1_2.log 2> err &
After this step you will need the WFK and SCR files produced in this run for the next runs. Move tgw1o_DS2_WFK to tgw1o_DS1_WFK and tgw1o_DS3_SCR to tgw1o_DS1_SCR.
The next sections are intended to show you how to find the converged parameters for a GW calculation. In principle, the following parameters might be used to decrease the CPU time and/or the memory requirements: optdriver = 3 ecuteps, nband and, for optdriver = 4, nband.
Before 2008, the advice was indeed to check independently what was the best value for each of these. However, with the evolution of memory/disk space, as well as the advent of new techniques to diminish the number of bands that is needed (see e.g. [Bruneval2008] and the input variable gwcomp), standard calculations nowadays only need the tuning of nband ecuteps, simultaneously for optdriver=3 and =4. Indeed, ecutwfn and can have the default value of ecut, while ecutsigx can have the default value of 4 * ecut for norm-conserving pseudopotentials, or pawecutdg for PAW calculations.
We begin by the convergence study on the only important parameter needed in the self- energy calculation (optdriver = 4): nband. This is because for these, we will not need a double dataset loop to check this convergence, and we will rely on the previously determined SCR file.
3 Convergence on the number of bands to calculate \Sigma_c¶
Let us check the convergence on the number of bands in the calculation of \Sigma_c. This convergence study is rather important, usually, BUT it can be done at the same time as the convergence study for the number of bands for the dielectric matrix.
The convergence on the number of bands to calculate the Self-Energy will be done by defining five datasets, with increasing nband:
ndtset 5 nband: 50 nband+ 50
In directory Work_gw1, copy the file ../tgw1_3.in, and modify the tgw1_x.files file as usual. Edit the tgw1_3.in file, and take the time to examine it. Then, issue:
abinit < tgw1_x.files > tgw1_3.log 2> err &
# Crystalline silicon # Calculation of the GW corrections # Convergence with respect to the number of bands ndtset 5 nband: 50 nband+ 50 # Calculation of the Self-Energy matrix elements (GW corrections) optdriver 4 getwfk 2 getscr 3 ecutsigx 8.0 ppmfrq 16.7 eV nkptgw 1 kptgw 0.000 0.000 0.000 bdgw 4 5 # Definition of the unit cell: fcc acell 3*10.217 # This is equivalent to 10.217 10.217 10.217 rprim 0.0 0.5 0.5 # FCC primitive vectors (to be scaled by acell) 0.5 0.0 0.5 0.5 0.5 0.0 # Definition of the atom types ntypat 1 # There is only one type of atom znucl 14 # The keyword "znucl" refers to the atomic number of the # possible type(s) of atom. The pseudopotential(s) # mentioned in the "files" file must correspond # to the type(s) of atom. Here, the only type is Silicon. # Definition of the atoms natom 2 # There are two atoms typat 1 1 # They both are of type 1, that is, Silicon. xred # Reduced coordinate of atoms 0.0 0.0 0.0 0.25 0.25 0.25 # Definition of the planewave basis set (at convergence 16 Rydberg 8 Hartree) ecut 8.0 # Maximal kinetic energy cut-off, in Hartree # Sampling of the BZ ngkpt 1 1 1 nshiftk 1 shiftk 0.0 0.0 0.0 istwfk *1 # This is mandatory in all the GW steps. #%%<BEGIN TEST_INFO> #%% [setup] #%% executable = abinit #%% test_chain = tgw1_2.in, tgw1_3.in, tgw1_4.in, tgw1_5.in #%% [files] #%% files_to_test = #%% tgw1_3.out, tolnlines= 70, tolabs= 9.000e-03, tolrel= 3.000e-02 #%% psp_files = 14si.pspnc #%% [paral_info] #%% max_nprocs = 4 #%% [extra_info] #%% authors = Unknown #%% keywords = GW #%% description = #%% Crystalline silicon #%% Calculation of the GW corrections #%%<END TEST_INFO>
Edit the output file. The number of bands used for the self-energy is mentioned in the fragments of output:
SIGMA fundamental parameters: PLASMON POLE MODEL number of plane-waves for SigmaX 283 number of plane-waves for SigmaC and W 169 number of plane-waves for wavefunctions 283 number of bands 50
Gathering the GW energies for each number of bands, one gets:
number of bands 50 4 5.915 -11.652 -17.103 4.738 0.786 -0.273 -12.212 -0.560 5.355 5 8.445 -9.700 -3.222 -6.448 0.798 -0.254 -9.676 0.024 8.470 number of bands 100 4 5.915 -11.652 -17.103 4.660 0.785 -0.274 -12.273 -0.620 5.295 5 8.445 -9.700 -3.222 -6.522 0.797 -0.255 -9.734 -0.034 8.411 number of bands 150 4 5.915 -11.652 -17.103 4.649 0.785 -0.274 -12.281 -0.629 5.286 5 8.445 -9.700 -3.222 -6.531 0.797 -0.255 -9.742 -0.042 8.403 number of bands 200 4 5.915 -11.652 -17.103 4.646 0.785 -0.274 -12.284 -0.632 5.284 5 8.445 -9.700 -3.222 -6.534 0.797 -0.255 -9.745 -0.044 8.401 number of bands 250 4 5.915 -11.652 -17.103 4.645 0.785 -0.274 -12.284 -0.632 5.283 5 8.445 -9.700 -3.222 -6.535 0.797 -0.255 -9.745 -0.045 8.400
So that nband = 100 can be considered converged within 0.01 eV.
With AbiPy , one can use the abicomp script provides to compare multiple SIGRES.nc files
Use the --expose
option to visualize of the QP gaps extracted from the different netcdf files:
$ abicomp.py sigres tgw1_3o_*_SIGRES.nc -e -sns Output of robot.get_dataframe(): nsppol qpgap nspinor nspden nband nkpt \ tgw1_3o_DS1_SIGRES.nc 1 3.114257 1 1 50 1 tgw1_3o_DS2_SIGRES.nc 1 3.116411 1 1 100 1 tgw1_3o_DS3_SIGRES.nc 1 3.116962 1 1 150 1 tgw1_3o_DS4_SIGRES.nc 1 3.117476 1 1 200 1 tgw1_3o_DS5_SIGRES.nc 1 3.117396 1 1 250 1 ecutwfn ecuteps ecutsigx scr_nband sigma_nband \ tgw1_3o_DS1_SIGRES.nc 7.563851 5.105599 7.563851 25 50 tgw1_3o_DS2_SIGRES.nc 7.563851 5.105599 7.563851 25 100 tgw1_3o_DS3_SIGRES.nc 7.563851 5.105599 7.563851 25 150 tgw1_3o_DS4_SIGRES.nc 7.563851 5.105599 7.563851 25 200 tgw1_3o_DS5_SIGRES.nc 7.563851 5.105599 7.563851 25 250 gwcalctyp scissor_ene tgw1_3o_DS1_SIGRES.nc 0 0.0 tgw1_3o_DS2_SIGRES.nc 0 0.0 tgw1_3o_DS3_SIGRES.nc 0 0.0 tgw1_3o_DS4_SIGRES.nc 0 0.0 tgw1_3o_DS5_SIGRES.nc 0 0.0
Invoking the script without options will open an ipython terminal to interact with the AbiPy robot.
Use the -nb
option to automatically generate a jupyter notebook that will open in your browser.
For further details about the API provided by SigRes Robots see the Sigres notebook
and the notebook with the GW lesson for GW calculations powered by AbiPy.
4 Convergence on the number of bands to calculate the screening (\epsilon^{-1})¶
Now, we come back to the calculation of the screening. Adequate convergence studies will couple the change of parameters for optdriver = 3 with a computation of the GW energy changes. One cannot rely on the convergence of the macroscopic dielectric constant to assess the convergence of the GW energies.
As a consequence, we will define a double loop over the datasets:
ndtset 10 udtset 5 2
The datasets 12,22,32,42 and 52, drive the computation of the GW energies:
# Calculation of the Self-Energy matrix elements (GW corrections) optdriver?2 4 getscr?2 -1 ecutsigx 8.0 nband?2 100
The datasets 11,21,31,41 and 51, drive the corresponding computation of the screening:
# Calculation of the screening (epsilon^-1 matrix) optdriver?1 3
In this latter series, we will have to vary the two different parameters ecuteps and nband.
Let us begin with nband. This convergence study is rather important. It can be done at the same time as the convergence study for the number of bands for the self-energy. Note that the number of bands used to calculate both the screening and the self-energy can be lowered by a large amount by resorting to the extrapolar technique (see the input variable gwcomp).
Second, we check the convergence on the number of bands in the calculation of the screening. This will be done by defining five datasets, with increasing nband:
nband11 25 nband21 50 nband31 100 nband41 150 nband51 200
In directory Work_gw1, copy the file ../tgw1_4.in, and modify the tgw1_x.files file as usual. Edit the tgw1_4.in file, and take the time to examine it.
Then, issue:
abinit < tgw1_x.files > tgw1_4.log 2> err &
# Crystalline silicon # Calculation of the GW corrections ndtset 10 udtset 5 2 # Convergence with respect to the number of bands to calculate epsilon^-1 nband11 25 nband21 50 nband31 100 nband41 150 nband51 200 # Calculation of the screening (epsilon^-1 matrix) optdriver?1 3 ecuteps 6.0 # Calculation of the Self-Energy matrix elements (GW corrections) optdriver?2 4 getscr?2 -1 ecutsigx 8.0 nband?2 100 nkptgw?2 1 kptgw?2 0.000 0.000 0.000 bdgw?2 4 5 getwfk 2 ppmfrq 16.7 eV # Definition of the unit cell: fcc acell 3*10.217 # This is equivalent to 10.217 10.217 10.217 rprim 0.0 0.5 0.5 # FCC primitive vectors (to be scaled by acell) 0.5 0.0 0.5 0.5 0.5 0.0 # Definition of the atom types ntypat 1 # There is only one type of atom znucl 14 # The keyword "znucl" refers to the atomic number of the # possible type(s) of atom. The pseudopotential(s) # mentioned in the "files" file must correspond # to the type(s) of atom. Here, the only type is Silicon. # Definition of the atoms natom 2 # There are two atoms typat 1 1 # They both are of type 1, that is, Silicon. xred # Reduced coordinate of atoms 0.0 0.0 0.0 0.25 0.25 0.25 # Definition of the planewave basis set (at convergence 16 Rydberg 8 Hartree) ecut 8.0 # Maximal kinetic energy cut-off, in Hartree # Sampling of the BZ ngkpt 1 1 1 nshiftk 1 shiftk 0.0 0.0 0.0 istwfk *1 # This is mandatory in all the GW steps. #%%<BEGIN TEST_INFO> #%% [setup] #%% executable = abinit #%% test_chain = tgw1_2.in, tgw1_3.in, tgw1_4.in, tgw1_5.in #%% [files] #%% files_to_test = #%% tgw1_4.out, tolnlines= 70, tolabs= 7.000e-02, tolrel= 3.000e-02 #%% psp_files = 14si.pspnc #%% [paral_info] #%% max_nprocs = 4 #%% [extra_info] #%% authors = Unknown #%% keywords = GW #%% description = #%% Crystalline silicon #%% Calculation of the GW corrections #%%<END TEST_INFO>
Edit the output file. The number of bands used for the wavefunctions in the computation of the screening is mentioned in the fragments of output:
EPSILON^-1 parameters (SCR file): dimension of the eps^-1 matrix on file 169 dimension of the eps^-1 matrix used 169 number of plane-waves for wavefunctions 283 number of bands 25
Gathering the macroscopic dielectric constant and GW energies for each number of bands, one gets:
number of bands 25 dielectric constant = 96.4962 dielectric constant without local fields = 140.5247 4 5.915 -11.652 -17.103 4.660 0.785 -0.274 -12.273 -0.620 5.295 5 8.445 -9.700 -3.222 -6.522 0.797 -0.255 -9.734 -0.034 8.411 number of bands 50 dielectric constant = 97.6590 dielectric constant without local fields = 140.5293 4 5.915 -11.652 -17.103 4.471 0.785 -0.274 -12.421 -0.768 5.147 5 8.445 -9.700 -3.222 -6.710 0.795 -0.257 -9.884 -0.184 8.261 number of bands 100 dielectric constant = 98.3494 dielectric constant without local fields = 140.5307 4 5.915 -11.652 -17.103 4.384 0.785 -0.273 -12.490 -0.838 5.078 5 8.445 -9.700 -3.222 -6.800 0.794 -0.259 -9.956 -0.255 8.190 number of bands 150 dielectric constant = 98.5074 dielectric constant without local fields = 140.5309 4 5.915 -11.652 -17.103 4.363 0.785 -0.274 -12.506 -0.854 5.062 5 8.445 -9.700 -3.222 -6.820 0.794 -0.259 -9.971 -0.271 8.174 number of bands 200 dielectric constant = 98.5227 dielectric constant without local fields = 140.5310 4 5.915 -11.652 -17.103 4.353 0.784 -0.275 -12.513 -0.860 5.055 5 8.445 -9.700 -3.222 -6.827 0.794 -0.259 -9.977 -0.277 8.168
So that the computation using 100 bands can be considered converged within 0.01 eV. Note that the value of nband that converges the dielectric matrix is usually of the same order of magnitude than the one that converges \Sigma_c.
5 Convergence on the dimension of the screening \epsilon^{-1} matrix¶
Then, we check the convergence on the number of plane waves in the calculation of the screening. This will be done by defining six datasets, with increasing ecuteps:
ecuteps:? 3.0 ecuteps+? 1.0
In directory Work_gw1, copy the file ../tgw1_5.in, and modify the tgw1_x.files file as usual. Edit the tgw1_5.in file, and take the time to examine it.
Then, issue:
abinit < tgw1_x.files > tgw1_5.log 2> err &
# Crystalline silicon # Calculation of the GW corrections ndtset 12 udtset 6 2 # Calculation of the screening (epsilon^-1 matrix) optdriver?1 3 nband?1 100 # Convergence with respect to the dimension of epsilon^-1 matrix ecuteps:? 3.0 ecuteps+? 1.0 # Calculation of the Self-Energy matrix elements (GW corrections) optdriver?2 4 getscr?2 -1 ecutsigx 8.0 nband?2 100 nkptgw 1 kptgw 0.000 0.000 0.000 bdgw 4 5 # GW calculation general parameters getwfk 2 ppmfrq 16.7 eV # Definition of the unit cell: fcc acell 3*10.217 # This is equivalent to 10.217 10.217 10.217 rprim 0.0 0.5 0.5 # FCC primitive vectors (to be scaled by acell) 0.5 0.0 0.5 0.5 0.5 0.0 # Definition of the atom types ntypat 1 # There is only one type of atom znucl 14 # The keyword "znucl" refers to the atomic number of the # possible type(s) of atom. The pseudopotential(s) # mentioned in the "files" file must correspond # to the type(s) of atom. Here, the only type is Silicon. # Definition of the atoms natom 2 # There are two atoms typat 1 1 # They both are of type 1, that is, Silicon. xred # Reduced coordinate of atoms 0.0 0.0 0.0 0.25 0.25 0.25 # Definition of the planewave basis set (at convergence 16 Rydberg 8 Hartree) ecut 8.0 # Maximal kinetic energy cut-off, in Hartree # Sampling of the BZ ngkpt 1 1 1 nshiftk 1 shiftk 0.0 0.0 0.0 istwfk *1 # This is mandatory in all the GW steps. #%%<BEGIN TEST_INFO> #%% [setup] #%% executable = abinit #%% test_chain = tgw1_2.in, tgw1_3.in, tgw1_4.in, tgw1_5.in #%% [files] #%% files_to_test = #%% tgw1_5.out, tolnlines= 70, tolabs= 7.000e-02, tolrel= 7.000e-02 #%% psp_files = 14si.pspnc #%% [paral_info] #%% max_nprocs = 4 #%% [extra_info] #%% authors = Unknown #%% keywords = GW #%% description = #%% Crystalline silicon #%% Calculation of the GW corrections #%%<END TEST_INFO>
Edit the output file. The number of bands used for the wavefunctions in the computation of the screening is mentioned in the fragments of output:
EPSILON^-1 parameters (SCR file): dimension of the eps^-1 matrix 59
Gathering the macroscopic dielectric constant and GW energies for each number of bands, one gets:
dimension of the eps^-1 matrix 59 dielectric constant = 99.2682 dielectric constant without local fields = 140.5307 4 5.915 -11.652 -17.103 4.560 0.788 -0.269 -12.354 -0.701 5.214 5 8.445 -9.700 -3.222 -6.792 0.795 -0.258 -9.949 -0.249 8.196 dimension of the eps^-1 matrix 113 dielectric constant = 98.4253 dielectric constant without local fields = 140.5307 4 5.915 -11.652 -17.103 4.427 0.785 -0.273 -12.456 -0.804 5.112 5 8.445 -9.700 -3.222 -6.799 0.794 -0.259 -9.955 -0.255 8.191 dimension of the eps^-1 matrix 137 dielectric constant = 98.4218 dielectric constant without local fields = 140.5307 4 5.915 -11.652 -17.103 4.403 0.785 -0.273 -12.475 -0.823 5.093 5 8.445 -9.700 -3.222 -6.798 0.794 -0.259 -9.954 -0.254 8.192 dimension of the eps^-1 matrix 169 dielectric constant = 98.3494 dielectric constant without local fields = 140.5307 4 5.915 -11.652 -17.103 4.384 0.785 -0.273 -12.490 -0.838 5.078 5 8.445 -9.700 -3.222 -6.800 0.794 -0.259 -9.956 -0.255 8.190 dimension of the eps^-1 matrix 259 dielectric constant = 98.3147 dielectric constant without local fields = 140.5307 4 5.915 -11.652 -17.103 4.373 0.785 -0.273 -12.499 -0.846 5.069 5 8.445 -9.700 -3.222 -6.800 0.794 -0.259 -9.955 -0.255 8.190 dimension of the eps^-1 matrix 283 dielectric constant = 98.3130 dielectric constant without local fields = 140.5307 4 5.915 -11.652 -17.103 4.371 0.785 -0.274 -12.499 -0.847 5.068 5 8.445 -9.700 -3.222 -6.800 0.794 -0.259 -9.955 -0.255 8.190
So that ecuteps = 6.0 (%npweps = 169) can be considered converged within 0.01 eV.
At this stage, we know that for the screening computation, we need ecuteps = 6.0 Ha and nband = 100.
Of course, until now, we have skipped the most difficult part of the convergence tests: the convergence in the number of k-points. It is as important to check the convergence on this parameter, than on the other ones. However, this might be very time consuming, since the CPU time scales as the square of the number of k-points (roughly), and the number of k-points can increase very rapidly from one possible grid to the next denser one. This is why we will leave this out of the present tutorial, and consider that we already know a sufficient k-point grid, for the last calculation.
As discussed in [Setten2017], the convergence study for k-points the number of bands and the cutoff energies can be decoupled in the sense that one can start from a reasonaby coarse k-mesh to find the converged values of nband, ecuteps, ecutsigx and then fix these values and look at the convergence with respect to the BZ mesh.
6 Calculation of the GW corrections for the band gap at \Gamma¶
Now we try to perform a GW calculation for a real problem: the calculation of the GW corrections for the direct band gap of bulk Silicon at the \Gamma point.
In directory Work_gw1, copy the file ../tgw1_6.in, and modify the tgw1_x.files file as usual. Then, edit the tgw1_6.in file, and, without examining it, comment the line
ngkpt 2 2 2 # Density of k points used for the automatic tests of the tutorial
and uncomment the line
#ngkpt 4 4 4 # Density of k points needed for a converged calculation
Then, issue:
abinit < tgw1_x.files > tgw1_6.log 2> err &
This job lasts about 3-4 minutes so it is worth to run it before the examination of the input file. Now, you can examine it.
# Crystalline silicon # Calculation of the GW correction to the direct band gap in Gamma # Dataset 1: ground state calculation # Dataset 2: calculation of the WFK file # Dataset 3: calculation of the screening (epsilon^-1 matrix for W) # Dataset 4: calculation of the Self-Energy matrix elements (GW corrections) ndtset 4 ngkpt 2 2 2 # Density of k points used for the automatic tests of the tutorial #ngkpt 4 4 4 # Density of k points needed for a converged calculation nshiftk 4 shiftk 0.0 0.0 0.0 # This grid contains the Gamma point, which is the point at which 0.0 0.5 0.5 # we will compute the (direct) band gap. There are 19 k points 0.5 0.0 0.5 # in the irreducible Brillouin zone, if ngkpt 4 4 4 is used. 0.5 0.5 0.0 istwfk *1 # For the GW computations, do not take advantage of the # specificities of k points to reduce the number of components of the # wavefunction. tolwfr 1.0d-10 # Dataset1: usual self-consistent ground-state calculation # Definition of the k-point grid shiftk1 0.5 0.5 0.5 # This grid is the most economical, but does not contain the Gamma point. 0.5 0.0 0.0 0.0 0.5 0.0 0.0 0.0 0.5 istwfk1 *0 # For the ground state, let the code use the time-reversal symmetry tolvrs1 1e-10 # Dataset2: calculation of WFK file # Definition of k-points iscf2 -2 # Non self-consistent calculation getden2 -1 # Read previous density file nband2 120 nbdbuf2 20 # Dataset3: Calculation of the screening (epsilon^-1 matrix) optdriver3 3 # Screening calculation getwfk3 -1 # Obtain WFK file from previous dataset nband3 25 # Bands to be used in the screening calculation ecuteps3 6.0 # Dimension of the screening matrix ppmfrq3 16.7 eV # Imaginary frequency where to calculate the screening # Dataset4: Calculation of the Self-Energy matrix elements (GW corrections) optdriver4 4 # Self-Energy calculation getwfk4 -2 # Obtain WFK file from dataset 1 getscr4 -1 # Obtain SCR file from previous dataset nband4 100 # Bands to be used in the Self-Energy calculation ecutsigx4 8.0 # Dimension of the G sum in Sigma_x nkptgw4 1 # number of k-point where to calculate the GW correction kptgw4 # k-points 0.000 0.000 0.000 # (Gamma) bdgw4 4 5 # calculate GW corrections for bands from 4 to 5 # Definition of the unit cell: fcc acell 3*10.217 # This is equivalent to 10.217 10.217 10.217 rprim 0.0 0.5 0.5 # FCC primitive vectors (to be scaled by acell) 0.5 0.0 0.5 0.5 0.5 0.0 # Definition of the atom types ntypat 1 # There is only one type of atom znucl 14 # The keyword "znucl" refers to the atomic number of the # possible type(s) of atom. The pseudopotential(s) # mentioned in the "files" file must correspond # to the type(s) of atom. Here, the only type is Silicon. # Definition of the atoms natom 2 # There are two atoms typat 1 1 # They both are of type 1, that is, Silicon. xred # Reduced coordinate of atoms 0.0 0.0 0.0 0.25 0.25 0.25 # Definition of the planewave basis set (at convergence 16 Rydberg 8 Hartree) ecut 8.0 # Maximal kinetic energy cut-off, in Hartree # Definition of the SCF procedure nstep 10 # Maximal number of SCF cycles diemac 12.0 # Although this is not mandatory, it is worth to # precondition the SCF cycle. The model dielectric # function used as the standard preconditioner # is described in the "dielng" input variable section. # Here, we follow the prescription for bulk silicon. #%%<BEGIN TEST_INFO> #%% [setup] #%% executable = abinit #%% [files] #%% files_to_test = #%% tgw1_6.out, tolnlines= 70, tolabs= 7.000e-02, tolrel= 3.000e-02, fld_options= -ridiculous #%% psp_files = 14si.pspnc #%% [paral_info] #%% max_nprocs = 4 #%% [extra_info] #%% authors = Unknown #%% keywords = GW #%% description = #%% Crystalline silicon #%% Calculation of the GW correction to the direct band gap in Gamma #%% Dataset 1: ground state calculation #%% Dataset 2: calculation of the WFK file #%% Dataset 3: calculation of the screening (epsilon^-1 matrix for W) #%% Dataset 4: calculation of the Self-Energy matrix elements (GW corrections) #%%<END TEST_INFO>
We need the usual part of the input file to perform a ground state calculation. This is done in dataset 1 and at the end we print out the density. We use a 4x4x4 FCC grid (so, 256 k-points in the full Brillouin Zone), shifted, because it is the most economical. It gives 10 k-points in the Irreducible part of the Brillouin Zone. However, this k-point grid does not contains the \Gamma point, and one cannot perform calculations of the self-energy corrections for other k-points than those present in the grid of k-points in the WFK file.
Then in dataset 2 we perform a non self-consistent calculation to calculate the KS structure in a set of 19 k-points in the Irreducible Brillouin Zone. This set of k-points is also derived from a 4x4x4 FCC grid, but a NON- SHIFTED one. It has the same density of points as the 10 k-point set, but the symmetries are not used in the most efficient way. However, this set contains the \Gamma point, which allows us to tackle the computation of the band gap at this point.
In dataset 3 we calculate the screening. The screening calculation is very time-consuming. So, we have decided to decrease a bit the parameters found in the previous convergence studies. Indeed, nband has been decreased from 100 to 25. This is a drastic change. The CPU time of this part is linear with respect to this parameter (or more exactly, with the number of conduction bands). Thus, the CPU time has been decreased by a factor of 4. Referring to our previous convergence study, we see that the absolute accuracy on the GW energies is now on the order of 0.2 eV only. This would be annoying for the absolute positioning of the band energy as required for band-offset or ionization potential of finite systems. However, as long as we are only interested in the gap energy that is fine enough.
Finally in dataset 4 we calculate the self-energy matrix element at \Gamma, using the previously determined parameters.
You should obtain the following results:
--- !SelfEnergy_ee kpoint : [ 0.000, 0.000, 0.000, ] spin : 1 KS_gap : 2.513 QP_gap : 3.140 Delta_QP_KS: 0.627 data: !SigmaeeData | Band E0 <VxcDFT> SigX SigC(E0) Z dSigC/dE Sig(E) E-E0 E 4 5.951 -11.271 -13.260 1.449 0.766 -0.305 -11.685 -0.414 5.537 5 8.464 -10.056 -5.572 -4.207 0.767 -0.304 -9.843 0.213 8.677 ...
So that the DFT energy gap in \Gamma is about 2.51 eV, while the GW correction is about 0.63 eV, so that the GW band gap found is 3.14 eV.
One can compare now what have been obtained to what one can get from the literature.
EXP 3.40 eV Landolt-Boernstein DFT (LDA) LDA 2.57 eV L. Hedin, Phys. Rev. 139, A796 (1965) LDA 2.57 eV M.S. Hybertsen and S. Louie, PRL 55, 1418 (1985) LDA (FLAPW) 2.55 eV N. Hamada, M. Hwang and A.J. Freeman, PRB 41, 3620 (1990) LDA (PAW) 2.53 eV B. Arnaud and M. Alouani, PRB 62, 4464 (2000) LDA 2.53 eV present work GW 3.27 eV M.S. Hybertsen and S. Louie, PRL 55, 1418 (1985) GW 3.35 eV M.S. Hybertsen and S. Louie, PRB 34, 5390 (1986) GW 3.30 eV R.W. Godby, M. Schlueter, L.J. Sham, PRB 37, 10159 (1988) GW (FLAPW) 3.30 eV N. Hamada, M. Hwang and A.J. Freeman, PRB 41, 3620 (1990) GW (FLAPW) 3.12 eV W. Ku and A.G. Eguiluz, PRL 89, 126401 (2002) GW 3.17 eV present work
The values are spread over an interval of 0.2 eV. They depend on the details of the calculation. In the case of pseudopotential calculations, they depend of course on the pseudopotential used. However, a GW result is hardly meaningful within 0.1 eV, in the present state of the art. But this goes also with the other source of inaccuracy, the choice of the pseudopotential, that can arrive up to even 0.2 eV. This can also be taken into account when choosing the level of accuracy for the convergence parameters in the GW calculation.
How to compute GW band structures¶
Finally, it is possible to calculate a full GW band plot of a system via interpolation. There are three possible techniques.
The first one is based on the use of Wannier functions to interpolate a few selected points in the IBZ obtained using the direct GW approach [Hamann2009]. You need to have the Wannier90 plug-in installed. See the directory tests/wannier90, test case 03, for an example of a file where a GW calculation is followed by the use of Wannier90.
# Self-consistent GW test including wannier90 interface for GW quasiparticles # This test is poorly converged (see GW and wannier90 tutorials). # Silicon structure acell 10.263 10.263 10.263 rprim 0.00 0.50 0.50 0.50 0.00 0.50 0.50 0.50 0.00 natom 2 xred 0.00 0.00 0.00 0.25 0.25 0.25 ntypat 1 typat 1 1 znucl 14.00 symmorphi 0 symsigma 0 # Parameters common to all runs ecut 6.00 ecutsigx 1.49923969 ecuteps 1.49923969 istwfk 8*1 ngkpt 4 4 4 nstep 100 nshiftk 1 shiftk 0.00 0.00 0.00 enunit 2 ndtset 7 gwpara 1 # Self-consistent run to get the density toldfe1 1.00d-6 # Non-self-consistent run to get all cg wavefunctions getden2 1 getwfk2 1 iscf2 -2 tolwfr2 1.0d-10 nband2 30 # Calculation of the dielectric matrix - iteration 1 optdriver3 3 gwcalctyp3 28 getwfk3 2 nband3 10 ecutwfn3 1.49923969E+00 Hartree awtr3 0 # Note : the default awtr 1 is better # Calculation of the model GW corrections - iteration 1 optdriver4 4 gwcalctyp4 28 getwfk4 2 getscr4 3 nband4 10 ecutwfn4 1.49923969E+00 Hartree icutcoul4 3 # old deprecated value of icutcoul, only used for legacy # Calculation of the dielectric matrix - iteration 2 optdriver5 3 gwcalctyp5 28 getwfk5 2 getqps5 4 nband5 10 ecutwfn5 1.49923969E+00 Hartree awtr5 0 # Note : the default awtr 1 is better # Calculation of the model GW corrections - iteration 2 optdriver6 4 gwcalctyp6 28 getwfk6 2 getqps6 4 getscr6 5 nband6 10 ecutwfn6 1.49923969E+00 Hartree icutcoul6 3 # old deprecated value of icutcoul, only used for legacy # Calculation of the quasiparticle Wannier functions getden7 1 getwfk7 2 # Must point to lda wavefunction file completely # consistent with _WFK file generated for GW getqps7 6 kptopt7 3 # Must use full-zone k mesh for wannier90 istwfk7 64*1 iscf7 -2 nstep7 0 # Irreducible-zone wavefunctions will be transformed # using symmetry operations to fill the full zone, # and must not be modified (for consistency with GW) tolwfr7 1.0d-10 # Dummy, but necessary nband7 10 # Must be consistent with nband in quasiparticle GW above prtwant7 3 # Flag for wannier90 interface with quaisparticles w90iniprj7 2 # Flag to use hydrogenic or gaussian orbital initial # projectors (to be specified in *.win file) w90prtunk7 0 # Flag for producing UNK files necessary for plotting # (suppressed here by 0 value) #Common to all model GW calculations rhoqpmix 0.5 nkptgw 8 kptgw 0.00000000E+00 0.00000000E+00 0.00000000E+00 2.50000000E-01 0.00000000E+00 0.00000000E+00 5.00000000E-01 0.00000000E+00 0.00000000E+00 2.50000000E-01 2.50000000E-01 0.00000000E+00 5.00000000E-01 2.50000000E-01 0.00000000E+00 -2.50000000E-01 2.50000000E-01 0.00000000E+00 5.00000000E-01 5.00000000E-01 0.00000000E+00 -2.50000000E-01 5.00000000E-01 2.50000000E-01 bdgw 1 8 # Only bands 1-8 are quasiparticle. LDA will be used for # bands 9 and 10 in the wannier90 calculation. 1 8 1 8 1 8 1 8 1 8 1 8 1 8 #%%<BEGIN TEST_INFO> #%% [setup] #%% executable = abinit #%% [files] #%% files_to_test = #%% t03.out, tolnlines = 27, tolabs = 7.0e-03, tolrel = 1.05e-00 #%% extra_inputs = t03o_DS7_w90.win #%% psp_files = 14si.pspnc #%% [paral_info] #%% max_nprocs = 1 #%% [extra_info] #%% authors = D. Hamann #%% keywords = GW #%% description = #%% Cannot be executed in parallel (mlwfovl_qp) #%% Si fcc, in primitive cell (2 atoms/cell) #%% Test of self-consistent model GW (2 iterations) following Faleev et al., #%% [PRL 93, 126406 (2004)] followed by construction of quasiparticle #%% maximally-localized wannier functions [Hamann & Vanderbilt, #%% arXiv:0810.3616v1 (cond-mat.mtrl-sci)]. Cutoffs are set for test- #%% acceptable speed, and the results are poorly converged. The input #%% file is sufficiently annotated to serve as a model. Note that well- #%% converged GW calculations are extremely time consuming, and in general #%% it is advisable to run the SCGW part separately on a parallel system, #%% and then run a separate serial job modeled on the last dataset, #%% substituting "irdwfk" and "irdqps" for "getwfk" and "getqps," with #%% appropriate links to the files produced in the serial run. Note that #%% the _DEN file from the first dataset is also needed as input, although #%% the discontinued "irdden" input variable is not needed or supported. #%% Note that acceptable names for the secondary input file needed by the #%% wannier90 library are wannier90.win, w90.win, abo_DSn_w90.win (ndtset #%% >0) and abo_w90.win (ndtset=0), where abo is the 4th line of the .files #%% file and n is the wannier dataset. #%% topics = Wannier #%%<END TEST_INFO>
The wannier interpolation is a very accurate method, can handle band crossings but it may require additional work to obtain well localized wannier functions. Another practical way follows from the fact that the QP energies, similarly to the KS eigenvalues, must fulfill the symmetry properties:
and
where \GG is a reciprocal lattice vector and S is a rotation of the point group of the crystal. Therefore it’s possible to employ the star-function interpolation by Shankland, Koelling and Wood [Euwema1969], [Koelling1986] in the improved version proposed by [Pickett1988] to fit the ab-initio results. This interpolation technique, by construction, passes through the initial points and satisfies the basic symmetry property of the band energies. It should be stressed, however, that this Fourier-based method can have problems in the presence of band crossings that may cause unphysical oscillations between the ab-initio points. To reduce this spurious effect, we suggest to interpolate the GW corrections instead of the GW energies. The corrections, indeed, are usually smoother in k-space and the resulting fit is more stable. A python example showing how to construct an energy-dependent scissors operator with AbiPy is available here
The third method uses that fact that the GW corrections are usually linear with the energy, for each group of bands. This is evident when reporting on a plot the GW correction with respect to the 0-order KS energy for each state. One can then simply correct the KS band structure at any point, by using a GW correction for the k-points where it has not been calculated explicitly, using a fit of the GW correction at a sparse set of points. A python example showing how to construct an energy-dependent scissors operator with AbiPy is available here.
Advanced features in the GW code¶
The user might switch to the second GW tutorial before coming back to the present section.
Calculations without using the Plasmon-Pole model¶
In order to circumvent the plasmon-pole model, the GW frequency convolution has to be calculated explicitly along the real axis. This is a tough job, since G and W have poles along the real axis. Therefore it is more convenient to use another path of integration along the imaginary axis plus the residues enclosed in the path.
Consequently, it is better to evaluate the screening for imaginary frequencies (to perform the integration) and also for real frequencies (to evaluate the contributions of the residues that may enter into the path of integration). The number of imaginary frequencies is set by the input variable nfreqim. The regular grid of real frequencies is determined by the input variables nfreqre, which sets the number of real frequencies, and freqremax, which indicates the maximum real frequency used.
The method is particularly suited to output the spectral function (contained in file out.sig). The grid of real frequencies used to calculate the spectral function is set by the number of frequencies (input variable nfreqsp) and by the maximum frequency calculated (input variable freqspmax).
Self-consistent calculations¶
The details in the implementation and the justification for the approximations retained can be found in [Bruneval2006]. The only added input variables are getqps and irdqps. These variables concerns the reading of the _QPS file, that contains the eigenvalues and the unitary transform matrices of a previous quasiparticle calculation. QPS stands for “QuasiParticle Structure”.
The only modified input variables for self-consistent calculations are gwcalctyp and bdgw. When the variable gwcalctyp is in between 0 and 9, The code calculates the quasiparticle energies only and does not output any QPS file (as in a standard GW run). When the variable gwcalctyp is in between 10 and 19, the code calculates the quasiparticle energies only and outputs them in a QPS file. When the variable gwcalctyp is in between 20 and 29, the code calculates the quasiparticle energies and wavefunctions and outputs them in a QPS file.
For a full self-consistency calculation, the quasiparticle wavefunctions are expanded in the basis set of the KS wavefunctions. The variable bdgw now indicates the size of all matrices to be calculated and diagonalized. The quasiparticle wavefunctions are consequently linear combinations of the KS wavefunctions in between the min and max values of bdgw.
A correct self-consistent calculation should consist of the following runs:
- 1) Self-consistent KS calculation: outputs a WFK file
- 2) Screening calculation (with KS inputs): outputs a SCR file
- 3) Sigma calculation (with KS inputs): outputs a QPS file
- 4) Screening calculation (with the WFK, and QPS file as an input): outputs a new SCR file
- 5) Sigma calculation (with the WFK, QPS and the new SCR files): outputs a new QPS file
- 6) Screening calculation (with the WFK, the new QPS file): outputs a newer SCR file
- 7) Sigma calculation (with the WFK, the newer QPS and SCR files): outputs a newer QPS
- ............ and so on, until the desired accuracy is reached
Note that for Hartree-Fock calculations a dummy screening is required for initialization reasons. Therefore, a correct HF calculations should look like
- 1) Self-consistent KS calculation: outputs a WFK file
- 2) Screening calculation using very low convergence parameters (with KS inputs): output a dummy SCR file
- 3) Sigma calculation (with KS inputs): outputs a QPS file
- 4) Sigma calculation (with the WFK and QPS files): outputs a new QPS file
- 5) Sigma calculation (with the WFK and the new QPS file): outputs a newer QPS file
- ............ and so on, until the desired accuracy is reached
In the case of a self-consistent calculation, the output is slightly more complex: For instance, iteration 2
k = 0.500 0.250 0.000 Band E_DFT <VxcDFT> E(N-1) <Hhartree> SigX SigC[E(N-1)] Z dSigC/dE Sig[E(N)] DeltaE E(N)_pert E(N)_diago 1 -3.422 -10.273 -3.761 6.847 -15.232 4.034 1.000 0.000 -11.198 -0.590 -4.351 -4.351 2 -0.574 -10.245 -0.850 9.666 -13.806 2.998 1.000 0.000 -10.807 -0.291 -1.141 -1.141 3 2.242 -9.606 2.513 11.841 -11.452 1.931 1.000 0.000 -9.521 -0.193 2.320 2.320 4 3.595 -10.267 4.151 13.866 -11.775 1.842 1.000 0.000 -9.933 -0.217 3.934 3.934 5 7.279 -8.804 9.916 16.078 -4.452 -1.592 1.000 0.000 -6.044 0.119 10.034 10.035 6 10.247 -9.143 13.462 19.395 -4.063 -1.775 1.000 0.000 -5.838 0.095 13.557 13.557 7 11.488 -9.704 15.159 21.197 -4.061 -1.863 1.000 0.000 -5.924 0.113 15.273 15.273 8 11.780 -9.180 15.225 20.958 -3.705 -1.893 1.000 0.000 -5.598 0.135 15.360 15.360 E^0_gap 3.684 E^GW_gap 5.764 DeltaE^GW_gap 2.080
The columns are
- Band: index of the band
- E_DFT: DFT eigenvalue
- VxcDFT: diagonal expectation value of the xc potential in between DFT bra and ket
- E(N-1): quasiparticle energy of the previous iteration (equal to DFT for the first iteration)
- Hhartree: diagonal expectation value of the Hartree Hamiltonian (equal to E_DFT - VxcDFT for the first iteration only)
- SigX: diagonal expectation value of the exchange self-energy
- SigC[E(N-1)]: diagonal expectation value of the correlation self-energy (evaluated for the energy of the preceeding iteration)
- Z: quasiparticle renormalization factor Z (taken equal to 1 in methods HF, SEX, COHSEX and model GW)
- dSigC/dE: Derivative of the correlation self-energy with respect to the energy
- Sig[E(N)]: Total self-energy for the new quasiparticle energy
- DeltaE: Energy difference with respect to the previous step
- E(N)_pert: QP energy as obtained by the usual perturbative method
- E(N)_diago: QP energy as obtained by the full diagonalization